13 #include "../fem/fem.hpp"
14 #include "../general/text.hpp"
18 #if defined(_MSC_VER) && (_MSC_VER < 1800)
20 #define copysign _copysign
32 input >> Order >> NumOfControlPoints;
33 knot.Load(input, NumOfControlPoints + Order + 1);
40 NumOfControlPoints = NCP;
41 knot.SetSize(NumOfControlPoints + Order + 1);
63 " Parent KnotVector order higher than child");
66 const int nOrder = Order +
t;
69 for (
int i = 0; i <= nOrder; i++)
71 (*newkv)[i] = knot(0);
73 for (
int i = nOrder + 1; i < newkv->
GetNCP(); i++)
75 (*newkv)[i] = knot(i - t);
77 for (
int i = 0; i <= nOrder; i++)
79 (*newkv)[newkv->
GetNCP() + i] = knot(knot.Size()-1);
89 newknots.
SetSize(NumOfElements);
91 for (
int i = 0; i < knot.Size()-1; i++)
93 if (knot(i) != knot(i+1))
95 newknots(j) = 0.5*(knot(i) + knot(i+1));
104 for (
int i = Order; i < NumOfControlPoints; i++)
106 if (knot(i) != knot(i+1))
115 double apb = knot(0) + knot(knot.Size()-1);
117 int ns = (NumOfControlPoints - Order)/2;
118 for (
int i = 1; i <= ns; i++)
120 double tmp = apb - knot(Order + i);
121 knot(Order + i) = apb - knot(NumOfControlPoints - i);
122 knot(NumOfControlPoints - i) = tmp;
128 os << Order <<
' ' << NumOfControlPoints <<
' ';
129 knot.Print(os, knot.Size());
137 double x, dx = 1.0/double (samples - 1);
139 for (
int i = 0; i <GetNE() ; i++)
141 for (
int j = 0; j <samples; j++)
146 CalcShape ( shape, i, x);
147 for (
int d = 0; d < Order+1; d++) { os<<
"\t"<<shape[d]; }
149 CalcDShape ( shape, i, x);
150 for (
int d = 0; d < Order+1; d++) { os<<
"\t"<<shape[d]; }
152 CalcD2Shape ( shape, i, x);
153 for (
int d = 0; d < Order+1; d++) { os<<
"\t"<<shape[d]; }
162 MFEM_ASSERT(Order <= MaxOrder, "Order > MaxOrder!
");
165 int ip = (i >= 0) ? (i + p) : (-1 - i + p);
166 double u = getKnotLocation((i >= 0) ? xi : 1. - xi, ip), saved, tmp;
167 double left[MaxOrder+1], right[MaxOrder+1];
170 for (int j = 1; j <= p; ++j)
172 left[j] = u - knot(ip+1-j);
173 right[j] = knot(ip+j) - u;
175 for (int r = 0; r < j; ++r)
177 tmp = shape(r)/(right[r+1] + left[j-r]);
178 shape(r) = saved + right[r+1]*tmp;
179 saved = left[j-r]*tmp;
185 // Routine from "The NURBS book
" - 2nd ed - Piegl and Tiller
186 void KnotVector::CalcDShape(Vector &grad, int i, double xi) const
188 int p = Order, rk, pk;
189 int ip = (i >= 0) ? (i + p) : (-1 - i + p);
190 double u = getKnotLocation((i >= 0) ? xi : 1. - xi, ip), temp, saved, d;
191 double ndu[MaxOrder+1][MaxOrder+1], left[MaxOrder+1], right[MaxOrder+1];
201 for (int j = 1; j <= p; j++)
203 left[j] = u - knot(ip-j+1);
204 right[j] = knot(ip+j) - u;
206 for (int r = 0; r < j; r++)
208 ndu[j][r] = right[r+1] + left[j-r];
209 temp = ndu[r][j-1]/ndu[j][r];
210 ndu[r][j] = saved + right[r+1]*temp;
211 saved = left[j-r]*temp;
216 for (int r = 0; r <= p; ++r)
223 d = ndu[rk][pk]/ndu[p][rk];
227 d -= ndu[r][pk]/ndu[p][r];
234 grad *= p*(knot(ip+1) - knot(ip));
238 grad *= p*(knot(ip) - knot(ip+1));
242 // Routine from "The NURBS book
" - 2nd ed - Piegl and Tiller
243 void KnotVector::CalcDnShape(Vector &gradn, int n, int i, double xi) const
245 int p = Order, rk, pk, j1, j2,r,j,k;
246 int ip = (i >= 0) ? (i + p) : (-1 - i + p);
247 double u = getKnotLocation((i >= 0) ? xi : 1. - xi, ip);
248 double temp, saved, d;
249 double a[2][MaxOrder+1],ndu[MaxOrder+1][MaxOrder+1], left[MaxOrder+1],
260 for (j = 1; j <= p; j++)
262 left[j] = u - knot(ip-j+1);
263 right[j] = knot(ip+j)- u;
266 for (r = 0; r < j; r++)
268 ndu[j][r] = right[r+1] + left[j-r];
269 temp = ndu[r][j-1]/ndu[j][r];
270 ndu[r][j] = saved + right[r+1]*temp;
271 saved = left[j-r]*temp;
276 for (r = 0; r <= p; r++)
281 for (k = 1; k <= n; k++)
288 a[s2][0] = a[s1][0]/ndu[pk+1][rk];
289 d = a[s2][0]*ndu[rk][pk];
310 for (j = j1; j <= j2; j++)
312 a[s2][j] = (a[s1][j] - a[s1][j-1])/ndu[pk+1][rk+j];
313 d += a[s2][j]*ndu[rk+j][pk];
318 a[s2][k] = - a[s1][k-1]/ndu[pk+1][r];
319 d += a[s2][j]*ndu[rk+j][pk];
330 u = (knot(ip+1) - knot(ip));
334 u = (knot(ip) - knot(ip+1));
338 for (k = 1; k <= n-1; k++) { temp *= (p-k)*u; }
340 for (j = 0; j <= p; j++) { gradn[j] *= temp; }
345 int KnotVector::findKnotSpan(double u) const
349 if (u == knot(NumOfControlPoints+Order))
351 mid = NumOfControlPoints;
356 high = NumOfControlPoints + 1;
357 mid = (low + high)/2;
358 while ( (u < knot(mid-1)) || (u > knot(mid)) )
368 mid = (low + high)/2;
374 void KnotVector::Difference(const KnotVector &kv, Vector &diff) const
376 if (Order != kv.GetOrder())
379 " Can not compare knot vectors with different orders!
");
382 int s = kv.Size() - Size();
385 kv.Difference(*this, diff);
393 for (int j = 0; j < kv.Size(); j++)
395 if (knot(i) == kv[j])
407 void NURBSPatch::init(int dim_)
414 ni = kv[0]->GetNCP();
415 nj = kv[1]->GetNCP();
418 data = new double[ni*nj*Dim];
421 for (int i = 0; i < ni*nj*Dim; i++)
427 else if (kv.Size() == 3)
429 ni = kv[0]->GetNCP();
430 nj = kv[1]->GetNCP();
431 nk = kv[2]->GetNCP();
433 data = new double[ni*nj*nk*Dim];
436 for (int i = 0; i < ni*nj*nk*Dim; i++)
448 NURBSPatch::NURBSPatch(const NURBSPatch &orig)
449 : ni(orig.ni), nj(orig.nj), nk(orig.nk), Dim(orig.Dim),
450 data(NULL), kv(orig.kv.Size()), sd(orig.sd), nd(orig.nd)
452 // Allocate and copy data:
453 const int data_size = Dim*ni*nj*((kv.Size() == 2) ? 1 : nk);
454 data = new double[data_size];
455 std::memcpy(data, orig.data, data_size*sizeof(double));
457 // Copy the knot vectors:
458 for (int i = 0; i < kv.Size(); i++)
460 kv[i] = new KnotVector(*orig.kv[i]);
464 NURBSPatch::NURBSPatch(std::istream &input)
466 int pdim, dim, size = 1;
469 input >> ws >> ident >> pdim; // knotvectors
471 for (int i = 0; i < pdim; i++)
473 kv[i] = new KnotVector(input);
474 size *= kv[i]->GetNCP();
477 input >> ws >> ident >> dim; // dimension
480 input >> ws >> ident; // controlpoints (homogeneous coordinates)
481 if (ident == "controlpoints
" || ident == "controlpoints_homogeneous
")
483 for (int j = 0, i = 0; i < size; i++)
484 for (int d = 0; d <= dim; d++, j++)
489 else // "controlpoints_cartesian
" (Cartesian coordinates with weight)
491 for (int j = 0, i = 0; i < size; i++)
493 for (int d = 0; d <= dim; d++)
497 for (int d = 0; d < dim; d++)
499 data[j+d] *= data[j+dim];
506 NURBSPatch::NURBSPatch(const KnotVector *kv0, const KnotVector *kv1, int dim_)
509 kv[0] = new KnotVector(*kv0);
510 kv[1] = new KnotVector(*kv1);
514 NURBSPatch::NURBSPatch(const KnotVector *kv0, const KnotVector *kv1,
515 const KnotVector *kv2, int dim_)
518 kv[0] = new KnotVector(*kv0);
519 kv[1] = new KnotVector(*kv1);
520 kv[2] = new KnotVector(*kv2);
524 NURBSPatch::NURBSPatch(Array<const KnotVector *> &kv_, int dim_)
526 kv.SetSize(kv_.Size());
527 for (int i = 0; i < kv.Size(); i++)
529 kv[i] = new KnotVector(*kv_[i]);
534 NURBSPatch::NURBSPatch(NURBSPatch *parent, int dir, int Order, int NCP)
536 kv.SetSize(parent->kv.Size());
537 for (int i = 0; i < kv.Size(); i++)
540 kv[i] = new KnotVector(*parent->kv[i]);
544 kv[i] = new KnotVector(Order, NCP);
549 void NURBSPatch::swap(NURBSPatch *np)
556 for (int i = 0; i < kv.Size(); i++)
558 if (kv[i]) { delete kv[i]; }
575 NURBSPatch::~NURBSPatch()
582 for (int i = 0; i < kv.Size(); i++)
584 if (kv[i]) { delete kv[i]; }
588 void NURBSPatch::Print(std::ostream &os) const
592 os << "knotvectors\n
" << kv.Size() << '\n';
593 for (int i = 0; i < kv.Size(); i++)
596 size *= kv[i]->GetNCP();
599 os << "\ndimension\n
" << Dim - 1
600 << "\n\ncontrolpoints\n
";
601 for (int j = 0, i = 0; i < size; i++)
604 for (int d = 1; d < Dim; d++)
606 os << ' ' << data[j++];
612 int NURBSPatch::SetLoopDirection(int dir)
633 " Direction error in 2D patch, dir =
" << dir << '\n';
663 " Direction error in 3D patch, dir =
" << dir << '\n';
671 void NURBSPatch::UniformRefinement()
674 for (int dir = 0; dir < kv.Size(); dir++)
676 kv[dir]->UniformRefinement(newknots);
677 KnotInsert(dir, newknots);
681 void NURBSPatch::KnotInsert(Array<KnotVector *> &newkv)
683 for (int dir = 0; dir < kv.Size(); dir++)
685 KnotInsert(dir, *newkv[dir]);
689 void NURBSPatch::KnotInsert(int dir, const KnotVector &newkv)
691 if (dir >= kv.Size() || dir < 0)
696 int t = newkv.GetOrder() - kv[dir]->GetOrder();
700 DegreeElevate(dir, t);
708 GetKV(dir)->Difference(newkv, diff);
711 KnotInsert(dir, diff);
715 void NURBSPatch::KnotInsert(Array<Vector *> &newkv)
717 for (int dir = 0; dir < kv.Size(); dir++)
719 KnotInsert(dir, *newkv[dir]);
723 // Routine from "The NURBS book
" - 2nd ed - Piegl and Tiller
724 void NURBSPatch::KnotInsert(int dir, const Vector &knot)
726 if (knot.Size() == 0 ) { return; }
728 if (dir >= kv.Size() || dir < 0)
733 NURBSPatch &oldp = *this;
734 KnotVector &oldkv = *kv[dir];
736 NURBSPatch *newpatch = new NURBSPatch(this, dir, oldkv.GetOrder(),
737 oldkv.GetNCP() + knot.Size());
738 NURBSPatch &newp = *newpatch;
739 KnotVector &newkv = *newp.GetKV(dir);
741 int size = oldp.SetLoopDirection(dir);
742 if (size != newp.SetLoopDirection(dir))
747 int rr = knot.Size() - 1;
748 int a = oldkv.findKnotSpan(knot(0)) - 1;
749 int b = oldkv.findKnotSpan(knot(rr)) - 1;
750 int pl = oldkv.GetOrder();
751 int ml = oldkv.GetNCP();
753 for (int j = 0; j <= a; j++)
757 for (int j = b+pl; j <= ml+pl; j++)
759 newkv[j+rr+1] = oldkv[j];
761 for (int k = 0; k <= (a-pl); k++)
763 for (int ll = 0; ll < size; ll++)
765 newp(k,ll) = oldp(k,ll);
768 for (int k = (b-1); k < ml; k++)
770 for (int ll = 0; ll < size; ll++)
772 newp(k+rr+1,ll) = oldp(k,ll);
779 for (int j = rr; j >= 0; j--)
781 while ( (knot(j) <= oldkv[i]) && (i > a) )
784 for (int ll = 0; ll < size; ll++)
786 newp(k-pl-1,ll) = oldp(i-pl-1,ll);
793 for (int ll = 0; ll < size; ll++)
795 newp(k-pl-1,ll) = newp(k-pl,ll);
798 for (int l = 1; l <= pl; l++)
801 double alfa = newkv[k+l] - knot(j);
802 if (fabs(alfa) == 0.0)
804 for (int ll = 0; ll < size; ll++)
806 newp(ind-1,ll) = newp(ind,ll);
811 alfa = alfa/(newkv[k+l] - oldkv[i-pl+l]);
812 for (int ll = 0; ll < size; ll++)
814 newp(ind-1,ll) = alfa*newp(ind-1,ll) + (1.0-alfa)*newp(ind,ll);
828 void NURBSPatch::DegreeElevate(int t)
830 for (int dir = 0; dir < kv.Size(); dir++)
832 DegreeElevate(dir, t);
836 // Routine from "The NURBS book
" - 2nd ed - Piegl and Tiller
837 void NURBSPatch::DegreeElevate(int dir, int t)
839 if (dir >= kv.Size() || dir < 0)
844 int i, j, k, kj, mpi, mul, mh, kind, cind, first, last;
845 int r, a, b, oldr, save, s, tr, lbz, rbz, l;
846 double inv, ua, ub, numer, alf, den, bet, gam;
848 NURBSPatch &oldp = *this;
849 KnotVector &oldkv = *kv[dir];
851 NURBSPatch *newpatch = new NURBSPatch(this, dir, oldkv.GetOrder() + t,
852 oldkv.GetNCP() + oldkv.GetNE()*t);
853 NURBSPatch &newp = *newpatch;
854 KnotVector &newkv = *newp.GetKV(dir);
856 int size = oldp.SetLoopDirection(dir);
857 if (size != newp.SetLoopDirection(dir))
862 int p = oldkv.GetOrder();
863 int n = oldkv.GetNCP()-1;
865 DenseMatrix bezalfs (p+t+1, p+1);
866 DenseMatrix bpts (p+1, size);
867 DenseMatrix ebpts (p+t+1, size);
868 DenseMatrix nextbpts(p-1, size);
876 Array2D<int> binom(ph+1, ph+1);
877 for (i = 0; i <= ph; i++)
879 binom(i,0) = binom(i,i) = 1;
880 for (j = 1; j < i; j++)
882 binom(i,j) = binom(i-1,j) + binom(i-1,j-1);
889 for (i = 1; i <= ph2; i++)
891 inv = 1.0/binom(ph,i);
893 for (j = max(0,i-t); j <= mpi; j++)
895 bezalfs(i,j) = inv*binom(p,j)*binom(t,i-j);
900 for (i = ph2+1; i < ph; i++)
903 for (j = max(0,i-t); j <= mpi; j++)
905 bezalfs(i,j) = bezalfs(ph-i,p-j);
916 for (l = 0; l < size; l++)
918 newp(0,l) = oldp(0,l);
920 for (i = 0; i <= ph; i++)
925 for (i = 0; i <= p; i++)
927 for (l = 0; l < size; l++)
929 bpts(i,l) = oldp(i,l);
936 while (b < m && oldkv[b] == oldkv[b+1]) { b++; }
944 if (oldr > 0) { lbz = (oldr+2)/2; }
947 if (r > 0) { rbz = ph-(r+1)/2; }
953 for (k = p ; k > mul; k--)
955 alphas[k-mul-1] = numer/(oldkv[a+k]-ua);
958 for (j = 1; j <= r; j++)
962 for (k = p; k >= s; k--)
964 for (l = 0; l < size; l++)
965 bpts(k,l) = (alphas[k-s]*bpts(k,l) +
966 (1.0-alphas[k-s])*bpts(k-1,l));
968 for (l = 0; l < size; l++)
970 nextbpts(save,l) = bpts(p,l);
975 for (i = lbz; i <= ph; i++)
977 for (l = 0; l < size; l++)
982 for (j = max(0,i-t); j <= mpi; j++)
984 for (l = 0; l < size; l++)
986 ebpts(i,l) += bezalfs(i,j)*bpts(j,l);
996 bet = (ub-newkv[kind-1])/den;
998 for (tr = 1; tr < oldr; tr++)
1007 alf = (ub-newkv[i])/(ua-newkv[i]);
1008 for (l = 0; l < size; l++)
1010 newp(i,l) = alf*newp(i,l)-(1.0-alf)*newp(i-1,l);
1015 if ((j-tr) <= (kind-ph+oldr))
1017 gam = (ub-newkv[j-tr])/den;
1018 for (l = 0; l < size; l++)
1020 ebpts(kj,l) = gam*ebpts(kj,l) + (1.0-gam)*ebpts(kj+1,l);
1025 for (l = 0; l < size; l++)
1027 ebpts(kj,l) = bet*ebpts(kj,l) + (1.0-bet)*ebpts(kj+1,l);
1042 for (i = 0; i < (ph-oldr); i++)
1048 for (j = lbz; j <= rbz; j++)
1050 for (l = 0; l < size; l++)
1052 newp(cind,l) = ebpts(j,l);
1059 for (j = 0; j <r; j++)
1060 for (l = 0; l < size; l++)
1062 bpts(j,l) = nextbpts(j,l);
1065 for (j = r; j <= p; j++)
1066 for (l = 0; l < size; l++)
1068 bpts(j,l) = oldp(b-p+j,l);
1077 for (i = 0; i <= ph; i++)
1083 newkv.GetElements();
1088 void NURBSPatch::FlipDirection(int dir)
1090 int size = SetLoopDirection(dir);
1092 for (int id = 0; id < nd/2; id++)
1093 for (int i = 0; i < size; i++)
1095 Swap<double>((*this)(id,i), (*this)(nd-1-id,i));
1100 void NURBSPatch::SwapDirections(int dir1, int dir2)
1102 if (abs(dir1-dir2) == 2)
1105 " directions 0 and 2 are not supported!
");
1108 Array<const KnotVector *> nkv(kv);
1110 Swap<const KnotVector *>(nkv[dir1], nkv[dir2]);
1111 NURBSPatch *newpatch = new NURBSPatch(nkv, Dim);
1113 int size = SetLoopDirection(dir1);
1114 newpatch->SetLoopDirection(dir2);
1116 for (int id = 0; id < nd; id++)
1117 for (int i = 0; i < size; i++)
1119 (*newpatch)(id,i) = (*this)(id,i);
1125 void NURBSPatch::Get3DRotationMatrix(double n[], double angle, double r,
1129 double l2 = n[0]*n[0] + n[1]*n[1] + n[2]*n[2];
1130 double l = sqrt(l2);
1132 if (fabs(angle) == M_PI_2)
1134 s = r*copysign(1., angle);
1138 else if (fabs(angle) == M_PI)
1153 T(0,0) = (n[0]*n[0] + (n[1]*n[1] + n[2]*n[2])*c)/l2;
1154 T(0,1) = -(n[0]*n[1]*c1)/l2 - (n[2]*s)/l;
1155 T(0,2) = -(n[0]*n[2]*c1)/l2 + (n[1]*s)/l;
1156 T(1,0) = -(n[0]*n[1]*c1)/l2 + (n[2]*s)/l;
1157 T(1,1) = (n[1]*n[1] + (n[0]*n[0] + n[2]*n[2])*c)/l2;
1158 T(1,2) = -(n[1]*n[2]*c1)/l2 - (n[0]*s)/l;
1159 T(2,0) = -(n[0]*n[2]*c1)/l2 - (n[1]*s)/l;
1160 T(2,1) = -(n[1]*n[2]*c1)/l2 + (n[0]*s)/l;
1161 T(2,2) = (n[2]*n[2] + (n[0]*n[0] + n[1]*n[1])*c)/l2;
1164 void NURBSPatch::Rotate3D(double n[], double angle)
1172 Vector x(3), y(NULL, 3);
1174 Get3DRotationMatrix(n, angle, 1., T);
1177 for (int i = 0; i < kv.Size(); i++)
1179 size *= kv[i]->GetNCP();
1182 for (int i = 0; i < size; i++)
1184 y.SetData(data + i*Dim);
1190 int NURBSPatch::MakeUniformDegree(int degree)
1196 for (int dir = 0; dir < kv.Size(); dir++)
1198 maxd = std::max(maxd, kv[dir]->GetOrder());
1202 for (int dir = 0; dir < kv.Size(); dir++)
1204 if (maxd > kv[dir]->GetOrder())
1206 DegreeElevate(dir, maxd - kv[dir]->GetOrder());
1213 NURBSPatch *Interpolate(NURBSPatch &p1, NURBSPatch &p2)
1215 if (p1.kv.Size() != p2.kv.Size() || p1.Dim != p2.Dim)
1220 int size = 1, dim = p1.Dim;
1221 Array<const KnotVector *> kv(p1.kv.Size() + 1);
1223 for (int i = 0; i < p1.kv.Size(); i++)
1225 if (p1.kv[i]->GetOrder() < p2.kv[i]->GetOrder())
1227 p1.KnotInsert(i, *p2.kv[i]);
1228 p2.KnotInsert(i, *p1.kv[i]);
1232 p2.KnotInsert(i, *p1.kv[i]);
1233 p1.KnotInsert(i, *p2.kv[i]);
1236 size *= kv[i]->GetNCP();
1239 KnotVector &nkv = *(new KnotVector(1, 2));
1240 nkv[0] = nkv[1] = 0.0;
1241 nkv[2] = nkv[3] = 1.0;
1245 NURBSPatch *patch = new NURBSPatch(kv, dim);
1248 for (int i = 0; i < size; i++)
1250 for (int d = 0; d < dim; d++)
1252 patch->data[i*dim+d] = p1.data[i*dim+d];
1253 patch->data[(i+size)*dim+d] = p2.data[i*dim+d];
1260 NURBSPatch *Revolve3D(NURBSPatch &patch, double n[], double ang, int times)
1268 Array<const KnotVector *> nkv(patch.kv.Size() + 1);
1270 for (int i = 0; i < patch.kv.Size(); i++)
1272 nkv[i] = patch.kv[i];
1273 size *= nkv[i]->GetNCP();
1276 KnotVector &lkv = *(new KnotVector(2, ns));
1278 lkv[0] = lkv[1] = lkv[2] = 0.0;
1279 for (int i = 1; i < times; i++)
1281 lkv[2*i+1] = lkv[2*i+2] = i;
1283 lkv[ns] = lkv[ns+1] = lkv[ns+2] = times;
1285 NURBSPatch *newpatch = new NURBSPatch(nkv, 4);
1288 DenseMatrix T(3), T2(3);
1289 Vector u(NULL, 3), v(NULL, 3);
1291 NURBSPatch::Get3DRotationMatrix(n, ang, 1., T);
1292 double c = cos(ang/2);
1293 NURBSPatch::Get3DRotationMatrix(n, ang/2, 1./c, T2);
1296 double *op = patch.data, *np;
1297 for (int i = 0; i < size; i++)
1299 np = newpatch->data + 4*i;
1300 for (int j = 0; j < 4; j++)
1304 for (int j = 0; j < times; j++)
1307 v.SetData(np += 4*size);
1310 v.SetData(np += 4*size);
1321 NURBSExtension::NURBSExtension(const NURBSExtension &orig)
1322 : mOrder(orig.mOrder), mOrders(orig.mOrders),
1323 NumOfKnotVectors(orig.NumOfKnotVectors),
1324 NumOfVertices(orig.NumOfVertices),
1325 NumOfElements(orig.NumOfElements),
1326 NumOfBdrElements(orig.NumOfBdrElements),
1327 NumOfDofs(orig.NumOfDofs),
1328 NumOfActiveVertices(orig.NumOfActiveVertices),
1329 NumOfActiveElems(orig.NumOfActiveElems),
1330 NumOfActiveBdrElems(orig.NumOfActiveBdrElems),
1331 NumOfActiveDofs(orig.NumOfActiveDofs),
1332 activeVert(orig.activeVert),
1333 activeElem(orig.activeElem),
1334 activeBdrElem(orig.activeBdrElem),
1335 activeDof(orig.activeDof),
1336 patchTopo(new Mesh(*orig.patchTopo)),
1338 edge_to_knot(orig.edge_to_knot),
1339 knotVectors(orig.knotVectors.Size()), // knotVectors are copied in the body
1340 weights(orig.weights),
1341 d_to_d(orig.d_to_d),
1342 master(orig.master),
1344 v_meshOffsets(orig.v_meshOffsets),
1345 e_meshOffsets(orig.e_meshOffsets),
1346 f_meshOffsets(orig.f_meshOffsets),
1347 p_meshOffsets(orig.p_meshOffsets),
1348 v_spaceOffsets(orig.v_spaceOffsets),
1349 e_spaceOffsets(orig.e_spaceOffsets),
1350 f_spaceOffsets(orig.f_spaceOffsets),
1351 p_spaceOffsets(orig.p_spaceOffsets),
1352 el_dof(orig.el_dof ? new Table(*orig.el_dof) : NULL),
1353 bel_dof(orig.bel_dof ? new Table(*orig.bel_dof) : NULL),
1354 el_to_patch(orig.el_to_patch),
1355 bel_to_patch(orig.bel_to_patch),
1356 el_to_IJK(orig.el_to_IJK),
1357 bel_to_IJK(orig.bel_to_IJK),
1358 patches(orig.patches.Size()) // patches are copied in the body
1360 // Copy the knot vectors:
1361 for (int i = 0; i < knotVectors.Size(); i++)
1363 knotVectors[i] = new KnotVector(*orig.knotVectors[i]);
1366 // Copy the patches:
1367 for (int p = 0; p < patches.Size(); p++)
1369 patches[p] = new NURBSPatch(*orig.patches[p]);
1373 NURBSExtension::NURBSExtension(std::istream &input)
1376 patchTopo = new Mesh;
1377 patchTopo->LoadPatchTopo(input, edge_to_knot);
1381 // CheckBdrPatches();
1383 skip_comment_lines(input, '#');
1385 // Read knotvectors or patches
1387 input >> ws >> ident; // 'knotvectors' or 'patches'
1388 if (ident == "knotvectors
")
1390 input >> NumOfKnotVectors;
1391 knotVectors.SetSize(NumOfKnotVectors);
1392 for (int i = 0; i < NumOfKnotVectors; i++)
1394 knotVectors[i] = new KnotVector(input);
1397 else if (ident == "patches
")
1399 patches.SetSize(GetNP());
1400 for (int p = 0; p < patches.Size(); p++)
1402 skip_comment_lines(input, '#');
1403 patches[p] = new NURBSPatch(input);
1406 NumOfKnotVectors = 0;
1407 for (int i = 0; i < patchTopo->GetNEdges(); i++)
1408 if (NumOfKnotVectors < KnotInd(i))
1410 NumOfKnotVectors = KnotInd(i);
1413 knotVectors.SetSize(NumOfKnotVectors);
1416 Array<int> edges, oedge;
1417 for (int p = 0; p < patches.Size(); p++)
1419 patchTopo->GetElementEdges(p, edges, oedge);
1420 if (Dimension() == 2)
1422 if (knotVectors[KnotInd(edges[0])] == NULL)
1424 knotVectors[KnotInd(edges[0])] =
1425 new KnotVector(*patches[p]->GetKV(0));
1427 if (knotVectors[KnotInd(edges[1])] == NULL)
1429 knotVectors[KnotInd(edges[1])] =
1430 new KnotVector(*patches[p]->GetKV(1));
1435 if (knotVectors[KnotInd(edges[0])] == NULL)
1437 knotVectors[KnotInd(edges[0])] =
1438 new KnotVector(*patches[p]->GetKV(0));
1440 if (knotVectors[KnotInd(edges[3])] == NULL)
1442 knotVectors[KnotInd(edges[3])] =
1443 new KnotVector(*patches[p]->GetKV(1));
1445 if (knotVectors[KnotInd(edges[8])] == NULL)
1447 knotVectors[KnotInd(edges[8])] =
1448 new KnotVector(*patches[p]->GetKV(2));
1455 MFEM_ABORT("invalid section:
" << ident);
1458 SetOrdersFromKnotVectors();
1463 // NumOfVertices, NumOfElements, NumOfBdrElements, NumOfDofs
1465 skip_comment_lines(input, '#');
1467 // Check for a list of mesh elements
1468 if (patches.Size() == 0)
1470 input >> ws >> ident;
1472 if (patches.Size() == 0 && ident == "mesh_elements
")
1474 input >> NumOfActiveElems;
1475 activeElem.SetSize(GetGNE());
1478 for (int i = 0; i < NumOfActiveElems; i++)
1481 activeElem[glob_elem] = true;
1484 skip_comment_lines(input, '#');
1485 input >> ws >> ident;
1489 NumOfActiveElems = NumOfElements;
1490 activeElem.SetSize(NumOfElements);
1494 GenerateActiveVertices();
1496 GenerateElementDofTable();
1497 GenerateActiveBdrElems();
1498 GenerateBdrElementDofTable();
1501 if (ident == "periodic
")
1506 skip_comment_lines(input, '#');
1507 input >> ws >> ident;
1510 if (patches.Size() == 0)
1513 if (ident == "weights
")
1515 weights.Load(input, GetNDof());
1517 else // e.g. ident = "unitweights
" or "autoweights
"
1519 weights.SetSize(GetNDof());
1525 ConnectBoundaries();
1528 NURBSExtension::NURBSExtension(NURBSExtension *parent, int newOrder)
1530 patchTopo = parent->patchTopo;
1533 parent->edge_to_knot.Copy(edge_to_knot);
1535 NumOfKnotVectors = parent->GetNKV();
1536 knotVectors.SetSize(NumOfKnotVectors);
1537 const Array<int> &pOrders = parent->GetOrders();
1538 for (int i = 0; i < NumOfKnotVectors; i++)
1540 if (newOrder > pOrders[i])
1543 parent->GetKnotVector(i)->DegreeElevate(newOrder - pOrders[i]);
1547 knotVectors[i] = new KnotVector(*parent->GetKnotVector(i));
1551 // copy some data from parent
1552 NumOfElements = parent->NumOfElements;
1553 NumOfBdrElements = parent->NumOfBdrElements;
1555 SetOrdersFromKnotVectors();
1557 GenerateOffsets(); // dof offsets will be different from parent
1559 NumOfActiveVertices = parent->NumOfActiveVertices;
1560 NumOfActiveElems = parent->NumOfActiveElems;
1561 NumOfActiveBdrElems = parent->NumOfActiveBdrElems;
1562 parent->activeVert.Copy(activeVert);
1564 parent->activeElem.Copy(activeElem);
1565 parent->activeBdrElem.Copy(activeBdrElem);
1567 GenerateElementDofTable();
1568 GenerateBdrElementDofTable();
1570 weights.SetSize(GetNDof());
1574 parent->master.Copy(master);
1575 parent->slave.Copy(slave);
1576 ConnectBoundaries();
1579 NURBSExtension::NURBSExtension(NURBSExtension *parent,
1580 const Array<int> &newOrders)
1582 newOrders.Copy(mOrders);
1583 SetOrderFromOrders();
1585 patchTopo = parent->patchTopo;
1588 parent->edge_to_knot.Copy(edge_to_knot);
1590 NumOfKnotVectors = parent->GetNKV();
1591 MFEM_VERIFY(mOrders.Size() == NumOfKnotVectors, "invalid newOrders array
");
1592 knotVectors.SetSize(NumOfKnotVectors);
1593 const Array<int> &pOrders = parent->GetOrders();
1595 for (int i = 0; i < NumOfKnotVectors; i++)
1597 if (mOrders[i] > pOrders[i])
1600 parent->GetKnotVector(i)->DegreeElevate(mOrders[i] - pOrders[i]);
1604 knotVectors[i] = new KnotVector(*parent->GetKnotVector(i));
1608 // copy some data from parent
1609 NumOfElements = parent->NumOfElements;
1610 NumOfBdrElements = parent->NumOfBdrElements;
1612 GenerateOffsets(); // dof offsets will be different from parent
1614 NumOfActiveVertices = parent->NumOfActiveVertices;
1615 NumOfActiveElems = parent->NumOfActiveElems;
1616 NumOfActiveBdrElems = parent->NumOfActiveBdrElems;
1617 parent->activeVert.Copy(activeVert);
1619 parent->activeElem.Copy(activeElem);
1620 parent->activeBdrElem.Copy(activeBdrElem);
1622 GenerateElementDofTable();
1623 GenerateBdrElementDofTable();
1625 weights.SetSize(GetNDof());
1628 parent->master.Copy(master);
1629 parent->slave.Copy(slave);
1630 ConnectBoundaries();
1633 NURBSExtension::NURBSExtension(Mesh *mesh_array[], int num_pieces)
1635 NURBSExtension *parent = mesh_array[0]->NURBSext;
1637 if (!parent->own_topo)
1640 " parent does not own the patch topology!
");
1642 patchTopo = parent->patchTopo;
1644 parent->own_topo = 0;
1646 parent->edge_to_knot.Copy(edge_to_knot);
1648 parent->GetOrders().Copy(mOrders);
1649 mOrder = parent->GetOrder();
1651 NumOfKnotVectors = parent->GetNKV();
1652 knotVectors.SetSize(NumOfKnotVectors);
1653 for (int i = 0; i < NumOfKnotVectors; i++)
1655 knotVectors[i] = new KnotVector(*parent->GetKnotVector(i));
1662 // assuming the meshes define a partitioning of all the elements
1663 NumOfActiveElems = NumOfElements;
1664 activeElem.SetSize(NumOfElements);
1667 GenerateActiveVertices();
1669 GenerateElementDofTable();
1670 GenerateActiveBdrElems();
1671 GenerateBdrElementDofTable();
1673 weights.SetSize(GetNDof());
1674 MergeWeights(mesh_array, num_pieces);
1677 NURBSExtension::~NURBSExtension()
1679 if (patches.Size() == 0)
1685 for (int i = 0; i < knotVectors.Size(); i++)
1687 delete knotVectors[i];
1690 for (int i = 0; i < patches.Size(); i++)
1701 void NURBSExtension::Print(std::ostream &os) const
1703 patchTopo->PrintTopo(os, edge_to_knot);
1704 if (patches.Size() == 0)
1706 os << "\nknotvectors\n
" << NumOfKnotVectors << '\n';
1707 for (int i = 0; i < NumOfKnotVectors; i++)
1709 knotVectors[i]->Print(os);
1712 if (NumOfActiveElems < NumOfElements)
1714 os << "\nmesh_elements\n
" << NumOfActiveElems << '\n';
1715 for (int i = 0; i < NumOfElements; i++)
1722 os << "\nweights\n
";
1723 weights.Print(os, 1);
1727 os << "\npatches\n
";
1728 for (int p = 0; p < patches.Size(); p++)
1730 os << "\n# patch
" << p << "\n\n
";
1731 patches[p]->Print(os);
1736 void NURBSExtension::PrintCharacteristics(std::ostream &os) const
1739 "NURBS
Mesh entity sizes:\n
"
1740 "Dimension =
" << Dimension() << "\n
"
1742 Array<int> unique_orders(mOrders);
1743 unique_orders.Sort();
1744 unique_orders.Unique();
1745 unique_orders.Print(os, unique_orders.Size());
1747 "NumOfKnotVectors =
" << GetNKV() << "\n
"
1748 "NumOfPatches =
" << GetNP() << "\n
"
1749 "NumOfBdrPatches =
" << GetNBP() << "\n
"
1750 "NumOfVertices =
" << GetGNV() << "\n
"
1751 "NumOfElements =
" << GetGNE() << "\n
"
1752 "NumOfBdrElements =
" << GetGNBE() << "\n
"
1753 "NumOfDofs =
" << GetNTotalDof() << "\n
"
1754 "NumOfActiveVertices =
" << GetNV() << "\n
"
1755 "NumOfActiveElems =
" << GetNE() << "\n
"
1756 "NumOfActiveBdrElems =
" << GetNBE() << "\n
"
1757 "NumOfActiveDofs =
" << GetNDof() << '\n';
1758 for (int i = 0; i < NumOfKnotVectors; i++)
1760 os << ' ' << i + 1 << ")
";
1761 knotVectors[i]->Print(os);
1766 void NURBSExtension::PrintFunctions(const char *basename, int samples) const
1769 for (int i = 0; i < NumOfKnotVectors; i++)
1771 std::ostringstream filename;
1772 filename << basename<<"_
"<<i<<".dat
";
1773 os.open(filename.str().c_str());
1774 knotVectors[i]->PrintFunctions(os,samples);
1779 void NURBSExtension::InitDofMap()
1786 void NURBSExtension::ConnectBoundaries(Array<int> &bnds0, Array<int> &bnds1)
1790 ConnectBoundaries();
1793 void NURBSExtension::ConnectBoundaries()
1795 if (master.Size() != slave.Size())
1799 if (master.Size() == 0 ) {
return; }
1802 d_to_d.SetSize(NumOfDofs);
1803 for (
int i = 0; i < NumOfDofs; i++) { d_to_d[i] = i; }
1805 for (
int i = 0; i < master.Size(); i++)
1807 if (Dimension() == 2)
1809 ConnectBoundaries2D(master[i], slave[i]);
1813 ConnectBoundaries3D(master[i], slave[i]);
1818 if (el_dof) {
delete el_dof; }
1819 if (bel_dof) {
delete bel_dof; }
1820 GenerateElementDofTable();
1821 GenerateBdrElementDofTable();
1826 int idx0 = -1, idx1 = -1;
1827 for (
int b = 0;
b < GetNBP();
b++)
1829 if (bnd0 == patchTopo->GetBdrAttribute(
b)) { idx0 =
b; }
1830 if (bnd1 == patchTopo->GetBdrAttribute(
b)) { idx1 =
b; }
1832 MFEM_VERIFY(idx0 != -1,
"Bdr 0 not found");
1833 MFEM_VERIFY(idx1 != -1,
"Bdr 1 not found");
1838 int okv0[1],okv1[1];
1845 int nks0 = kv0[0]->
GetNKS();
1848 bool compatible =
true;
1849 if (p2g0.
nx() != p2g1.
nx()) { compatible =
false; }
1850 if (kv0[0]->GetNKS() != kv1[0]->
GetNKS()) { compatible =
false; }
1851 if (kv0[0]->GetOrder() != kv1[0]->
GetOrder()) { compatible =
false; }
1858 mfem_error(
"NURBS boundaries not compatible");
1862 for (
int i = 0; i < nks0; i++)
1864 if (kv0[0]->isElement(i))
1866 if (!kv1[0]->isElement(i)) {
mfem_error(
"isElement does not match"); }
1867 for (
int ii = 0; ii <= kv0[0]->
GetOrder(); ii++)
1869 int ii0 = (okv0[0] >= 0) ? (i+ii) : (nx-i-ii);
1870 int ii1 = (okv1[0] >= 0) ? (i+ii) : (nx-i-ii);
1872 d_to_d[p2g0(ii0)] = d_to_d[p2g1(ii1)];
1882 for (
int i = 0; i < d_to_d.Size(); i++)
1888 for (
int i = 0; i < tmp.
Size(); i++)
1890 if (tmp[i] == 1) { tmp[i] = cnt++; }
1894 for (
int i = 0; i < d_to_d.Size(); i++)
1896 d_to_d[i] = tmp[d_to_d[i]];
1902 int idx0 = -1, idx1 = -1;
1903 for (
int b = 0;
b < GetNBP();
b++)
1905 if (bnd0 == patchTopo->GetBdrAttribute(
b)) { idx0 =
b; }
1906 if (bnd1 == patchTopo->GetBdrAttribute(
b)) { idx1 =
b; }
1908 MFEM_VERIFY(idx0 != -1,
"Bdr 0 not found");
1909 MFEM_VERIFY(idx1 != -1,
"Bdr 1 not found");
1914 int okv0[2],okv1[2];
1923 int nks0 = kv0[0]->
GetNKS();
1924 int nks1 = kv0[1]->
GetNKS();
1927 bool compatible =
true;
1928 if (p2g0.
nx() != p2g1.
nx()) { compatible =
false; }
1929 if (p2g0.
ny() != p2g1.
ny()) { compatible =
false; }
1931 if (kv0[0]->GetNKS() != kv1[0]->
GetNKS()) { compatible =
false; }
1932 if (kv0[1]->GetNKS() != kv1[0]->
GetNKS()) { compatible =
false; }
1934 if (kv0[0]->GetOrder() != kv1[0]->
GetOrder()) { compatible =
false; }
1935 if (kv0[1]->GetOrder() != kv1[1]->
GetOrder()) { compatible =
false; }
1947 mfem_error(
"NURBS boundaries not compatible");
1951 for (
int j = 0; j < nks1; j++)
1953 if (kv0[1]->isElement(j))
1955 if (!kv1[1]->isElement(j)) {
mfem_error(
"isElement does not match #1"); }
1956 for (
int i = 0; i < nks0; i++)
1958 if (kv0[0]->isElement(i))
1960 if (!kv1[0]->isElement(i)) {
mfem_error(
"isElement does not match #0"); }
1961 for (
int jj = 0; jj <= kv0[1]->
GetOrder(); jj++)
1963 int jj0 = (okv0[1] >= 0) ? (j+jj) : (ny-j-jj);
1964 int jj1 = (okv1[1] >= 0) ? (j+jj) : (ny-j-jj);
1966 for (
int ii = 0; ii <= kv0[0]->
GetOrder(); ii++)
1968 int ii0 = (okv0[0] >= 0) ? (i+ii) : (nx-i-ii);
1969 int ii1 = (okv1[0] >= 0) ? (i+ii) : (nx-i-ii);
1971 d_to_d[p2g0(ii0,jj0)] = d_to_d[p2g1(ii1,jj1)];
1983 for (
int i = 0; i < d_to_d.Size(); i++)
1989 for (
int i = 0; i < tmp.
Size(); i++)
1991 if (tmp[i] == 1) { tmp[i] = cnt++; }
1995 for (
int i = 0; i < d_to_d.Size(); i++)
1997 d_to_d[i] = tmp[d_to_d[i]];
2003 int vert[8], nv, g_el, nx, ny, nz,
dim = Dimension();
2009 activeVert.SetSize(GetGNV());
2011 for (
int p = 0;
p < GetNP();
p++)
2017 nz = (dim == 3) ? p2g.
nz() : 1;
2019 for (
int k = 0; k < nz; k++)
2021 for (
int j = 0; j < ny; j++)
2023 for (
int i = 0; i < nx; i++)
2025 if (activeElem[g_el])
2029 vert[0] = p2g(i, j );
2030 vert[1] = p2g(i+1,j );
2031 vert[2] = p2g(i+1,j+1);
2032 vert[3] = p2g(i, j+1);
2037 vert[0] = p2g(i, j, k);
2038 vert[1] = p2g(i+1,j, k);
2039 vert[2] = p2g(i+1,j+1,k);
2040 vert[3] = p2g(i, j+1,k);
2042 vert[4] = p2g(i, j, k+1);
2043 vert[5] = p2g(i+1,j, k+1);
2044 vert[6] = p2g(i+1,j+1,k+1);
2045 vert[7] = p2g(i, j+1,k+1);
2049 for (
int v = 0; v < nv; v++)
2051 activeVert[vert[v]] = 1;
2060 NumOfActiveVertices = 0;
2061 for (
int i = 0; i < GetGNV(); i++)
2062 if (activeVert[i] == 1)
2064 activeVert[i] = NumOfActiveVertices++;
2070 int dim = Dimension();
2073 activeBdrElem.SetSize(GetGNBE());
2074 if (GetGNE() == GetNE())
2076 activeBdrElem =
true;
2077 NumOfActiveBdrElems = GetGNBE();
2080 activeBdrElem =
false;
2081 NumOfActiveBdrElems = 0;
2094 for (
int i = 0; i < num_pieces; i++)
2100 for (
int lel = 0; lel < lext->
GetNE(); lel++)
2102 int gel = lelem_elem[lel];
2104 int nd = el_dof->RowSize(gel);
2105 int *gdofs = el_dof->GetRow(gel);
2107 for (
int j = 0; j < nd; j++)
2109 weights(gdofs[j]) = lext->
weights(ldofs[j]);
2122 for (
int i = 0; i < num_pieces; i++)
2129 for (
int lel = 0; lel < lext->
GetNE(); lel++)
2145 for (
int p = 0;
p < GetNP();
p++)
2147 patchTopo->GetElementEdges(
p, edges, oedge);
2149 for (
int i = 0; i < edges.
Size(); i++)
2151 edges[i] = edge_to_knot[edges[i]];
2154 edges[i] = -1 - edges[i];
2158 if ((Dimension() == 2 &&
2159 (edges[0] != -1 - edges[2] || edges[1] != -1 - edges[3])) ||
2161 (Dimension() == 3 &&
2162 (edges[0] != edges[2] || edges[0] != edges[4] ||
2163 edges[0] != edges[6] || edges[1] != edges[3] ||
2164 edges[1] != edges[5] || edges[1] != edges[7] ||
2165 edges[8] != edges[9] || edges[8] != edges[10] ||
2166 edges[8] != edges[11])))
2168 mfem::err <<
"NURBSExtension::CheckPatch (patch = " <<
p
2169 <<
")\n Inconsistent edge-to-knot mapping!\n";
2173 if ((Dimension() == 2 &&
2174 (edges[0] < 0 || edges[1] < 0)) ||
2176 (Dimension() == 3 &&
2177 (edges[0] < 0 || edges[3] < 0 || edges[8] < 0)))
2179 mfem::err <<
"NURBSExtension::CheckPatch (patch = " <<
p
2180 <<
") : Bad orientation!\n";
2191 for (
int p = 0;
p < GetNBP();
p++)
2193 patchTopo->GetBdrElementEdges(
p, edges, oedge);
2195 for (
int i = 0; i < edges.
Size(); i++)
2197 edges[i] = edge_to_knot[edges[i]];
2200 edges[i] = -1 - edges[i];
2204 if ((Dimension() == 2 && (edges[0] < 0)) ||
2205 (Dimension() == 3 && (edges[0] < 0 || edges[1] < 0)))
2207 mfem::err <<
"NURBSExtension::CheckBdrPatch (boundary patch = "
2208 <<
p <<
") : Bad orientation!\n";
2219 patchTopo->GetElementEdges(p, edges, orient);
2221 if (Dimension() == 2)
2223 kv[0] = KnotVec(edges[0]);
2224 kv[1] = KnotVec(edges[1]);
2228 kv[0] = KnotVec(edges[0]);
2229 kv[1] = KnotVec(edges[3]);
2230 kv[2] = KnotVec(edges[8]);
2240 patchTopo->GetElementEdges(p, edges, orient);
2242 if (Dimension() == 2)
2244 kv[0] = KnotVec(edges[0]);
2245 kv[1] = KnotVec(edges[1]);
2249 kv[0] = KnotVec(edges[0]);
2250 kv[1] = KnotVec(edges[3]);
2251 kv[2] = KnotVec(edges[8]);
2261 patchTopo->GetBdrElementEdges(p, edges, orient);
2263 if (Dimension() == 2)
2265 kv[0] = KnotVec(edges[0]);
2269 kv[0] = KnotVec(edges[0]);
2270 kv[1] = KnotVec(edges[1]);
2281 patchTopo->GetBdrElementEdges(p, edges, orient);
2283 if (Dimension() == 2)
2285 kv[0] = KnotVec(edges[0]);
2289 kv[0] = KnotVec(edges[0]);
2290 kv[1] = KnotVec(edges[1]);
2296 MFEM_VERIFY(mOrders.Size() > 0,
"");
2297 mOrder = mOrders[0];
2298 for (
int i = 1; i < mOrders.Size(); i++)
2300 if (mOrders[i] != mOrder)
2310 mOrders.SetSize(NumOfKnotVectors);
2311 for (
int i = 0; i < NumOfKnotVectors; i++)
2313 mOrders[i] = knotVectors[i]->GetOrder();
2315 SetOrderFromOrders();
2320 int nv = patchTopo->GetNV();
2321 int ne = patchTopo->GetNEdges();
2322 int nf = patchTopo->GetNFaces();
2323 int np = patchTopo->GetNE();
2324 int meshCounter, spaceCounter,
dim = Dimension();
2330 e_meshOffsets.SetSize(ne);
2331 f_meshOffsets.SetSize(nf);
2332 p_meshOffsets.SetSize(np);
2334 v_spaceOffsets.SetSize(nv);
2335 e_spaceOffsets.SetSize(ne);
2336 f_spaceOffsets.SetSize(nf);
2337 p_spaceOffsets.SetSize(np);
2340 for (meshCounter = 0; meshCounter < nv; meshCounter++)
2342 v_meshOffsets[meshCounter] = meshCounter;
2343 v_spaceOffsets[meshCounter] = meshCounter;
2345 spaceCounter = meshCounter;
2348 for (
int e = 0; e < ne; e++)
2350 e_meshOffsets[e] = meshCounter;
2351 e_spaceOffsets[e] = spaceCounter;
2352 meshCounter += KnotVec(e)->GetNE() - 1;
2353 spaceCounter += KnotVec(e)->GetNCP() - 2;
2357 for (
int f = 0;
f < nf;
f++)
2359 f_meshOffsets[
f] = meshCounter;
2360 f_spaceOffsets[
f] = spaceCounter;
2362 patchTopo->GetFaceEdges(
f, edges, orient);
2365 (KnotVec(edges[0])->GetNE() - 1) *
2366 (KnotVec(edges[1])->GetNE() - 1);
2368 (KnotVec(edges[0])->GetNCP() - 2) *
2369 (KnotVec(edges[1])->GetNCP() - 2);
2373 for (
int p = 0;
p < np;
p++)
2375 p_meshOffsets[
p] = meshCounter;
2376 p_spaceOffsets[
p] = spaceCounter;
2378 patchTopo->GetElementEdges(
p, edges, orient);
2383 (KnotVec(edges[0])->GetNE() - 1) *
2384 (KnotVec(edges[1])->GetNE() - 1);
2386 (KnotVec(edges[0])->GetNCP() - 2) *
2387 (KnotVec(edges[1])->GetNCP() - 2);
2392 (KnotVec(edges[0])->GetNE() - 1) *
2393 (KnotVec(edges[3])->GetNE() - 1) *
2394 (KnotVec(edges[8])->GetNE() - 1);
2396 (KnotVec(edges[0])->GetNCP() - 2) *
2397 (KnotVec(edges[3])->GetNCP() - 2) *
2398 (KnotVec(edges[8])->GetNCP() - 2);
2401 NumOfVertices = meshCounter;
2402 NumOfDofs = spaceCounter;
2407 int dim = Dimension();
2411 for (
int p = 0;
p < GetNP();
p++)
2413 GetPatchKnotVectors(
p, kv);
2415 int ne = kv[0]->GetNE();
2416 for (
int d = 1; d <
dim; d++)
2418 ne *= kv[d]->GetNE();
2421 NumOfElements += ne;
2427 int dim = Dimension() - 1;
2430 NumOfBdrElements = 0;
2431 for (
int p = 0;
p < GetNBP();
p++)
2433 GetBdrPatchKnotVectors(
p, kv);
2435 int ne = kv[0]->GetNE();
2436 for (
int d = 1; d <
dim; d++)
2438 ne *= kv[d]->GetNE();
2441 NumOfBdrElements += ne;
2449 if (Dimension() == 2)
2451 Get2DElementTopo(elements);
2455 Get3DElementTopo(elements);
2467 for (
int p = 0;
p < GetNP();
p++)
2473 int patch_attr = patchTopo->GetAttribute(
p);
2475 for (
int j = 0; j < ny; j++)
2477 for (
int i = 0; i < nx; i++)
2481 ind[0] = activeVert[p2g(i, j )];
2482 ind[1] = activeVert[p2g(i+1,j )];
2483 ind[2] = activeVert[p2g(i+1,j+1)];
2484 ind[3] = activeVert[p2g(i, j+1)];
2503 for (
int p = 0;
p < GetNP();
p++)
2510 int patch_attr = patchTopo->GetAttribute(
p);
2512 for (
int k = 0; k < nz; k++)
2514 for (
int j = 0; j < ny; j++)
2516 for (
int i = 0; i < nx; i++)
2520 ind[0] = activeVert[p2g(i, j, k)];
2521 ind[1] = activeVert[p2g(i+1,j, k)];
2522 ind[2] = activeVert[p2g(i+1,j+1,k)];
2523 ind[3] = activeVert[p2g(i, j+1,k)];
2525 ind[4] = activeVert[p2g(i, j, k+1)];
2526 ind[5] = activeVert[p2g(i+1,j, k+1)];
2527 ind[6] = activeVert[p2g(i+1,j+1,k+1)];
2528 ind[7] = activeVert[p2g(i, j+1,k+1)];
2530 elements[el] =
new Hexahedron(ind, patch_attr);
2544 if (Dimension() == 2)
2546 Get2DBdrElementTopo(boundary);
2550 Get3DBdrElementTopo(boundary);
2562 for (
int b = 0;
b < GetNBP();
b++)
2567 int bdr_patch_attr = patchTopo->GetBdrAttribute(
b);
2569 for (
int i = 0; i < nx; i++)
2571 if (activeBdrElem[g_be])
2573 int i_ = (okv[0] >= 0) ? i : (nx - 1 - i);
2574 ind[0] = activeVert[p2g[i_ ]];
2575 ind[1] = activeVert[p2g[i_+1]];
2577 boundary[l_be] =
new Segment(ind, bdr_patch_attr);
2593 for (
int b = 0;
b < GetNBP();
b++)
2599 int bdr_patch_attr = patchTopo->GetBdrAttribute(
b);
2601 for (
int j = 0; j < ny; j++)
2603 int j_ = (okv[1] >= 0) ? j : (ny - 1 - j);
2604 for (
int i = 0; i < nx; i++)
2606 if (activeBdrElem[g_be])
2608 int i_ = (okv[0] >= 0) ? i : (nx - 1 - i);
2609 ind[0] = activeVert[p2g(i_, j_ )];
2610 ind[1] = activeVert[p2g(i_+1,j_ )];
2611 ind[2] = activeVert[p2g(i_+1,j_+1)];
2612 ind[3] = activeVert[p2g(i_, j_+1)];
2625 activeDof.SetSize(GetNTotalDof());
2628 if (Dimension() == 2)
2630 Generate2DElementDofTable();
2634 Generate3DElementDofTable();
2637 NumOfActiveDofs = 0;
2638 for (
int d = 0; d < GetNTotalDof(); d++)
2642 activeDof[d] = NumOfActiveDofs;
2645 int *dof = el_dof->GetJ();
2646 int ndof = el_dof->Size_of_connections();
2647 for (
int i = 0; i < ndof; i++)
2649 dof[i] = activeDof[dof[i]] - 1;
2661 el_to_patch.
SetSize(NumOfActiveElems);
2662 el_to_IJK.SetSize(NumOfActiveElems, 2);
2664 for (
int p = 0;
p < GetNP();
p++)
2669 const int ord0 = kv[0]->
GetOrder();
2670 const int ord1 = kv[1]->
GetOrder();
2671 for (
int j = 0; j < kv[1]->
GetNKS(); j++)
2673 if (kv[1]->isElement(j))
2675 for (
int i = 0; i < kv[0]->
GetNKS(); i++)
2677 if (kv[0]->isElement(i))
2682 for (
int jj = 0; jj <= ord1; jj++)
2684 for (
int ii = 0; ii <= ord0; ii++)
2686 conn.
to = DofMap(p2g(i+ii,j+jj));
2687 activeDof[conn.
to] = 1;
2688 el_dof_list.
Append(conn);
2691 el_to_patch[el] =
p;
2692 el_to_IJK(el,0) = i;
2693 el_to_IJK(el,1) = j;
2704 el_dof =
new Table(NumOfActiveElems, el_dof_list);
2715 el_to_patch.
SetSize(NumOfActiveElems);
2716 el_to_IJK.SetSize(NumOfActiveElems, 3);
2718 for (
int p = 0;
p < GetNP();
p++)
2723 const int ord0 = kv[0]->
GetOrder();
2724 const int ord1 = kv[1]->
GetOrder();
2725 const int ord2 = kv[2]->
GetOrder();
2726 for (
int k = 0; k < kv[2]->
GetNKS(); k++)
2728 if (kv[2]->isElement(k))
2730 for (
int j = 0; j < kv[1]->
GetNKS(); j++)
2732 if (kv[1]->isElement(j))
2734 for (
int i = 0; i < kv[0]->
GetNKS(); i++)
2736 if (kv[0]->isElement(i))
2741 for (
int kk = 0; kk <= ord2; kk++)
2743 for (
int jj = 0; jj <= ord1; jj++)
2745 for (
int ii = 0; ii <= ord0; ii++)
2747 conn.
to = DofMap(p2g(i+ii, j+jj, k+kk));
2748 activeDof[conn.
to] = 1;
2749 el_dof_list.
Append(conn);
2754 el_to_patch[el] =
p;
2755 el_to_IJK(el,0) = i;
2756 el_to_IJK(el,1) = j;
2757 el_to_IJK(el,2) = k;
2770 el_dof =
new Table(NumOfActiveElems, el_dof_list);
2775 if (Dimension() == 2)
2777 Generate2DBdrElementDofTable();
2781 Generate3DBdrElementDofTable();
2784 int *dof = bel_dof->GetJ();
2785 int ndof = bel_dof->Size_of_connections();
2786 for (
int i = 0; i < ndof; i++)
2788 dof[i] = activeDof[dof[i]] - 1;
2795 int lbe = 0, okv[1];
2800 bel_to_patch.
SetSize(NumOfActiveBdrElems);
2801 bel_to_IJK.SetSize(NumOfActiveBdrElems, 1);
2803 for (
int b = 0;
b < GetNBP();
b++)
2806 const int nx = p2g.
nx();
2809 const int nks0 = kv[0]->
GetNKS();
2810 const int ord0 = kv[0]->
GetOrder();
2811 for (
int i = 0; i < nks0; i++)
2813 if (kv[0]->isElement(i))
2815 if (activeBdrElem[gbe])
2818 for (
int ii = 0; ii <= ord0; ii++)
2820 conn.
to = DofMap(p2g[(okv[0] >= 0) ? (i+ii) : (nx-i-ii)]);
2821 bel_dof_list.
Append(conn);
2823 bel_to_patch[lbe] =
b;
2824 bel_to_IJK(lbe,0) = (okv[0] >= 0) ? i : (-1-i);
2832 bel_dof =
new Table(NumOfActiveBdrElems, bel_dof_list);
2838 int lbe = 0, okv[2];
2843 bel_to_patch.
SetSize(NumOfActiveBdrElems);
2844 bel_to_IJK.SetSize(NumOfActiveBdrElems, 2);
2846 for (
int b = 0;
b < GetNBP();
b++)
2849 const int nx = p2g.
nx();
2850 const int ny = p2g.
ny();
2853 const int nks0 = kv[0]->
GetNKS();
2854 const int ord0 = kv[0]->
GetOrder();
2855 const int nks1 = kv[1]->
GetNKS();
2856 const int ord1 = kv[1]->
GetOrder();
2857 for (
int j = 0; j < nks1; j++)
2859 if (kv[1]->isElement(j))
2861 for (
int i = 0; i < nks0; i++)
2863 if (kv[0]->isElement(i))
2865 if (activeBdrElem[gbe])
2868 for (
int jj = 0; jj <= ord1; jj++)
2870 const int jj_ = (okv[1] >= 0) ? (j+jj) : (ny-j-jj);
2871 for (
int ii = 0; ii <= ord0; ii++)
2873 const int ii_ = (okv[0] >= 0) ? (i+ii) : (nx-i-ii);
2874 conn.
to = DofMap(p2g(ii_, jj_));
2875 bel_dof_list.
Append(conn);
2878 bel_to_patch[lbe] =
b;
2879 bel_to_IJK(lbe,0) = (okv[0] >= 0) ? i : (-1-i);
2880 bel_to_IJK(lbe,1) = (okv[1] >= 0) ? j : (-1-j);
2890 bel_dof =
new Table(NumOfActiveBdrElems, bel_dof_list);
2896 for (
int gv = 0; gv < GetGNV(); gv++)
2897 if (activeVert[gv] >= 0)
2899 lvert_vert[activeVert[gv]] = gv;
2906 for (
int le = 0, ge = 0; ge < GetGNE(); ge++)
2909 lelem_elem[le++] = ge;
2921 NURBSFE->
SetIJK(el_to_IJK.GetRow(i));
2922 if (el_to_patch[i] != NURBSFE->
GetPatch())
2924 GetPatchKnotVectors(el_to_patch[i], NURBSFE->
KnotVectors());
2928 el_dof->GetRow(i, dofs);
2929 weights.GetSubVector(dofs, NURBSFE->
Weights());
2942 NURBSFE->
SetIJK(bel_to_IJK.GetRow(i));
2943 if (bel_to_patch[i] != NURBSFE->
GetPatch())
2945 GetBdrPatchKnotVectors(bel_to_patch[i], NURBSFE->
KnotVectors());
2946 NURBSFE->
SetPatch(bel_to_patch[i]);
2949 bel_dof->GetRow(i, dofs);
2950 weights.GetSubVector(dofs, NURBSFE->
Weights());
2960 if (patches.Size() == 0)
2962 GetPatchNets(Nodes, Dimension());
2968 if (patches.Size() == 0) {
return; }
2970 SetSolutionVector(Nodes, Dimension());
2976 if (patches.Size() == 0)
2978 mfem_error(
"NURBSExtension::SetKnotsFromPatches :"
2979 " No patches available!");
2984 for (
int p = 0;
p < patches.Size();
p++)
2986 GetPatchKnotVectors(
p, kv);
2988 for (
int i = 0; i < kv.
Size(); i++)
2990 *kv[i] = *patches[
p]->GetKV(i);
2994 SetOrdersFromKnotVectors();
3001 NumOfActiveElems = NumOfElements;
3002 activeElem.SetSize(NumOfElements);
3005 GenerateActiveVertices();
3007 GenerateElementDofTable();
3008 GenerateActiveBdrElems();
3009 GenerateBdrElementDofTable();
3011 ConnectBoundaries();
3023 const int vdim = fes->
GetVDim();
3025 for (
int p = 0;
p < GetNP();
p++)
3030 const int nx = kv[0]->GetNCP();
3031 const int ny = kv[1]->GetNCP();
3032 const int nz = (kv.
Size() == 2) ? 1 : kv[2]->GetNCP();
3033 for (
int k = 0; k < nz; k++)
3035 for (
int j = 0; j < ny; j++)
3037 for (
int i = 0; i < nx; i++)
3039 const int ll = (kv.
Size() == 2) ? p2g(i,j) : p2g(i,j,k);
3040 const int l = DofMap(ll);
3041 for (
int vd = 0; vd < vdim; vd++)
3059 const int vdim = fes->
GetVDim();
3061 for (
int p = 0;
p < GetNP();
p++)
3063 os <<
"\n# patch " <<
p <<
"\n\n";
3066 const int nx = kv[0]->GetNCP();
3067 const int ny = kv[1]->GetNCP();
3068 const int nz = (kv.Size() == 2) ? 1 : kv[2]->GetNCP();
3069 for (
int k = 0; k < nz; k++)
3071 for (
int j = 0; j < ny; j++)
3073 for (
int i = 0; i < nx; i++)
3075 const int ll = (kv.Size() == 2) ? p2g(i,j) : p2g(i,j,k);
3076 const int l = DofMap(ll);
3078 for (
int vd = 1; vd < vdim; vd++)
3091 for (
int p = 0;
p < patches.Size();
p++)
3093 for (
int dir = 0; dir < patches[
p]->GetNKV(); dir++)
3095 int oldd = patches[
p]->GetKV(dir)->GetOrder();
3096 int newd = std::min(oldd + rel_degree, degree);
3099 patches[
p]->DegreeElevate(dir, newd - oldd);
3107 for (
int p = 0;
p < patches.Size();
p++)
3109 patches[
p]->UniformRefinement();
3120 for (
int p = 0;
p < patches.Size();
p++)
3122 patchTopo->GetElementEdges(
p, edges, orient);
3126 pkv[0] = kv[KnotInd(edges[0])];
3127 pkv[1] = kv[KnotInd(edges[1])];
3131 pkv[0] = kv[KnotInd(edges[0])];
3132 pkv[1] = kv[KnotInd(edges[3])];
3133 pkv[2] = kv[KnotInd(edges[8])];
3136 patches[
p]->KnotInsert(pkv);
3147 for (
int p = 0;
p < patches.Size();
p++)
3149 patchTopo->GetElementEdges(
p, edges, orient);
3153 pkv[0] = kv[KnotInd(edges[0])];
3154 pkv[1] = kv[KnotInd(edges[1])];
3158 pkv[0] = kv[KnotInd(edges[0])];
3159 pkv[1] = kv[KnotInd(edges[3])];
3160 pkv[2] = kv[KnotInd(edges[8])];
3163 patches[
p]->KnotInsert(pkv);
3170 if (Dimension() == 2)
3172 Get2DPatchNets(coords, vdim);
3176 Get3DPatchNets(coords, vdim);
3185 patches.SetSize(GetNP());
3186 for (
int p = 0;
p < GetNP();
p++)
3192 for (
int j = 0; j < kv[1]->GetNCP(); j++)
3194 for (
int i = 0; i < kv[0]->GetNCP(); i++)
3196 const int l = DofMap(p2g(i,j));
3197 for (
int d = 0; d < vdim; d++)
3199 Patch(i,j,d) = coords(l*vdim + d)*weights(l);
3201 Patch(i,j,vdim) = weights(l);
3212 patches.SetSize(GetNP());
3213 for (
int p = 0;
p < GetNP();
p++)
3219 for (
int k = 0; k < kv[2]->GetNCP(); k++)
3221 for (
int j = 0; j < kv[1]->GetNCP(); j++)
3223 for (
int i = 0; i < kv[0]->GetNCP(); i++)
3225 const int l = DofMap(p2g(i,j,k));
3226 for (
int d = 0; d < vdim; d++)
3228 Patch(i,j,k,d) = coords(l*vdim + d)*weights(l);
3230 Patch(i,j,k,vdim) = weights(l);
3239 if (Dimension() == 2)
3241 Set2DSolutionVector(coords, vdim);
3245 Set3DSolutionVector(coords, vdim);
3254 weights.SetSize(GetNDof());
3255 for (
int p = 0;
p < GetNP();
p++)
3259 MFEM_ASSERT(vdim+1 == Patch.
GetNC(),
"");
3261 for (
int j = 0; j < kv[1]->GetNCP(); j++)
3263 for (
int i = 0; i < kv[0]->GetNCP(); i++)
3265 const int l = p2g(i,j);
3266 for (
int d = 0; d < vdim; d++)
3268 coords(l*vdim + d) = Patch(i,j,d)/Patch(i,j,vdim);
3270 weights(l) = Patch(i,j,vdim);
3282 weights.SetSize(GetNDof());
3283 for (
int p = 0;
p < GetNP();
p++)
3287 MFEM_ASSERT(vdim+1 == Patch.
GetNC(),
"");
3289 for (
int k = 0; k < kv[2]->GetNCP(); k++)
3291 for (
int j = 0; j < kv[1]->GetNCP(); j++)
3293 for (
int i = 0; i < kv[0]->GetNCP(); i++)
3295 const int l = p2g(i,j,k);
3296 for (
int d = 0; d < vdim; d++)
3298 coords(l*vdim + d) = Patch(i,j,k,d)/Patch(i,j,k,vdim);
3300 weights(l) = Patch(i,j,k,vdim);
3312 partitioning(orig.partitioning ? new int[orig.GetGNE()] : NULL),
3314 ldof_group(orig.ldof_group)
3319 std::memcpy(partitioning, orig.partitioning, orig.
GetGNE()*
sizeof(int));
3331 mfem_error(
"ParNURBSExtension::ParNURBSExtension :\n"
3332 " all elements in the parent must be active!");
3339 mfem_error(
"ParNURBSExtension::ParNURBSExtension :\n"
3340 " parent does not own the patch topology!");
3362 partitioning =
new int[
GetGNE()];
3363 for (
int i = 0; i <
GetGNE(); i++)
3365 partitioning[i] = part[i];
3367 SetActive(partitioning, active_bel);
3375 BuildGroups(partitioning, *serial_elem_dof);
3379 for (
int gel = 0, lel = 0; gel <
GetGNE(); gel++)
3385 int *gdofs = serial_elem_dof->
GetRow(gel);
3386 for (
int i = 0; i < ndofs; i++)
3397 : gtopo(par_parent->gtopo.GetComm())
3455 partitioning = NULL;
3457 MFEM_VERIFY(par_parent->partitioning,
3458 "parent ParNURBSExtension has no partitioning!");
3463 bool extract_weights =
false;
3468 SetActive(par_parent->partitioning, par_parent->
activeBdrElem);
3479 extract_weights =
true;
3482 Table *glob_elem_dof = GetGlobalElementDofTable();
3483 BuildGroups(par_parent->partitioning, *glob_elem_dof);
3484 if (extract_weights)
3491 for (
int gel = 0, lel = 0; gel <
GetGNE(); gel++)
3497 int *gdofs = glob_elem_dof->
GetRow(gel);
3498 for (
int i = 0; i < ndofs; i++)
3500 weights(ldofs[i]) = glob_weights(gdofs[i]);
3506 delete glob_elem_dof;
3509 Table *ParNURBSExtension::GetGlobalElementDofTable()
3513 return Get2DGlobalElementDofTable();
3517 return Get3DGlobalElementDofTable();
3521 Table *ParNURBSExtension::Get2DGlobalElementDofTable()
3524 const KnotVector *kv[2];
3526 Array<Connection> gel_dof_list;
3530 p2g.SetPatchDofMap(
p, kv);
3533 const int ord0 = kv[0]->GetOrder();
3534 const int ord1 = kv[1]->GetOrder();
3535 for (
int j = 0; j < kv[1]->GetNKS(); j++)
3537 if (kv[1]->isElement(j))
3539 for (
int i = 0; i < kv[0]->GetNKS(); i++)
3541 if (kv[0]->isElement(i))
3543 Connection conn(el,0);
3544 for (
int jj = 0; jj <= ord1; jj++)
3546 for (
int ii = 0; ii <= ord0; ii++)
3548 conn.to = p2g(i+ii,j+jj);
3549 gel_dof_list.Append(conn);
3559 return (
new Table(
GetGNE(), gel_dof_list));
3562 Table *ParNURBSExtension::Get3DGlobalElementDofTable()
3565 const KnotVector *kv[3];
3567 Array<Connection> gel_dof_list;
3571 p2g.SetPatchDofMap(
p, kv);
3574 const int ord0 = kv[0]->GetOrder();
3575 const int ord1 = kv[1]->GetOrder();
3576 const int ord2 = kv[2]->GetOrder();
3577 for (
int k = 0; k < kv[2]->GetNKS(); k++)
3579 if (kv[2]->isElement(k))
3581 for (
int j = 0; j < kv[1]->GetNKS(); j++)
3583 if (kv[1]->isElement(j))
3585 for (
int i = 0; i < kv[0]->GetNKS(); i++)
3587 if (kv[0]->isElement(i))
3589 Connection conn(el,0);
3590 for (
int kk = 0; kk <= ord2; kk++)
3592 for (
int jj = 0; jj <= ord1; jj++)
3594 for (
int ii = 0; ii <= ord0; ii++)
3596 conn.to = p2g(i+ii,j+jj,k+kk);
3597 gel_dof_list.Append(conn);
3610 return (
new Table(
GetGNE(), gel_dof_list));
3613 void ParNURBSExtension::SetActive(
const int *partitioning_,
3614 const Array<bool> &active_bel)
3620 for (
int i = 0; i <
GetGNE(); i++)
3621 if (partitioning_[i] == MyRank)
3629 for (
int i = 0; i <
GetGNBE(); i++)
3636 void ParNURBSExtension::BuildGroups(
const int *partitioning_,
3637 const Table &elem_dof)
3641 ListOfIntegerSets groups;
3646 for (
int i = 0; i < dof_proc.Size_of_connections(); i++)
3648 dof_proc.GetJ()[i] = partitioning_[dof_proc.GetJ()[i]];
3653 group.Recreate(1, &MyRank);
3654 groups.Insert(group);
3661 group.Recreate(dof_proc.RowSize(d), dof_proc.GetRow(d));
3669 #endif // MFEM_USE_MPI
3672 void NURBSPatchMap::GetPatchKnotVectors(
int p,
const KnotVector *kv[])
3678 kv[0] = Ext->
KnotVec(edges[0]);
3679 kv[1] = Ext->
KnotVec(edges[1]);
3685 kv[0] = Ext->
KnotVec(edges[0]);
3686 kv[1] = Ext->
KnotVec(edges[3]);
3687 kv[2] = Ext->
KnotVec(edges[8]);
3692 void NURBSPatchMap::GetBdrPatchKnotVectors(
int p,
const KnotVector *kv[],
3697 kv[0] = Ext->
KnotVec(edges[0], oedge[0], &okv[0]);
3702 kv[1] = Ext->
KnotVec(edges[1], oedge[1], &okv[1]);
3712 GetPatchKnotVectors(p, kv);
3714 I = kv[0]->
GetNE() - 1;
3715 J = kv[1]->
GetNE() - 1;
3717 for (
int i = 0; i < verts.
Size(); i++)
3722 for (
int i = 0; i < edges.
Size(); i++)
3729 K = kv[2]->
GetNE() - 1;
3731 for (
int i = 0; i < faces.
Size(); i++)
3742 GetPatchKnotVectors(p, kv);
3747 for (
int i = 0; i < verts.
Size(); i++)
3752 for (
int i = 0; i < edges.
Size(); i++)
3761 for (
int i = 0; i < faces.
Size(); i++)
3773 GetBdrPatchKnotVectors(p, kv, okv);
3775 I = kv[0]->
GetNE() - 1;
3777 for (
int i = 0; i < verts.
Size(); i++)
3788 J = kv[1]->
GetNE() - 1;
3790 for (
int i = 0; i < edges.
Size(); i++)
3801 GetBdrPatchKnotVectors(p, kv, okv);
3805 for (
int i = 0; i < verts.
Size(); i++)
3818 for (
int i = 0; i < edges.
Size(); i++)
Array< KnotVector * > knotVectors
Abstract class for all finite elements.
Array< int > f_meshOffsets
KnotVector * DegreeElevate(int t) const
void Create(ListOfIntegerSets &groups, int mpitag)
Set up the group topology given the list of sets of shared entities.
Array2D< int > bel_to_IJK
int Size() const
Return the logical size of the array.
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
void SetCoordsFromPatches(Vector &Nodes)
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
void ConnectBoundaries2D(int bnd0, int bnd1)
KnotVector * KnotVec(int edge)
Class for grid function - Vector with associated FE space.
int DofToVDof(int dof, int vd, int ndofs=-1) const
void GetBdrElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of bdr element i.
void DegreeElevate(int dir, int t)
NURBSPatch * Interpolate(NURBSPatch &p1, NURBSPatch &p2)
void LoadBE(int i, const FiniteElement *BE) const
void Generate2DElementDofTable()
void MergeGridFunctions(GridFunction *gf_array[], int num_pieces, GridFunction &merged)
Array< const KnotVector * > & KnotVectors() const
void Generate3DElementDofTable()
void GetBdrElementFace(int i, int *f, int *o) const
Return the index and the orientation of the face of bdr element i. (3D)
void KnotInsert(int dir, const KnotVector &knot)
void MergeWeights(Mesh *mesh_array[], int num_pieces)
void SetSize(int s)
Resize the vector to size s.
void Get3DBdrElementTopo(Array< Element * > &boundary) const
void Copy(Array ©) const
Create a copy of the internal array to the provided copy.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
void Rotate3D(double normal[], double angle)
Array< int > e_meshOffsets
void GetElements()
Count the number of elements.
void SetOrdersFromKnotVectors()
An arbitrary order and dimension NURBS element.
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
void Generate3DBdrElementDofTable()
Array< int > e_spaceOffsets
void GenerateBdrElementDofTable()
const KnotVector * GetKnotVector(int i) const
void GetElementLocalToGlobal(Array< int > &lelem_elem)
Data type quadrilateral element.
Array< int > edge_to_knot
void GenerateActiveVertices()
void GetVertexLocalToGlobal(Array< int > &lvert_vert)
void Set2DSolutionVector(Vector &Nodes, int vdim)
void skip_comment_lines(std::istream &is, const char comment_char)
Check if the stream starts with comment_char. If so skip it.
void DeleteAll()
Delete the whole array.
void SwapDirections(int dir1, int dir2)
friend class NURBSPatchMap
KnotVector()
Create KnotVector.
static const int MaxOrder
Array< int > v_spaceOffsets
void UniformRefinement(Vector &newknots) const
void Generate2DBdrElementDofTable()
Data type hexahedron element.
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
void ConnectBoundaries3D(int bnd0, int bnd1)
ParNURBSExtension(const ParNURBSExtension &orig)
void SetIJK(const int *IJK) const
double f(const Vector &xvec)
void Set3DSolutionVector(Vector &Nodes, int vdim)
void SetPatch(int p) const
void GetPatchKnotVectors(int p, Array< KnotVector * > &kv)
void GenerateActiveBdrElems()
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Mesh * GetMesh() const
Returns the mesh.
void PrintSolution(const GridFunction &sol, std::ostream &out) const
void Get3DPatchNets(const Vector &Nodes, int vdim)
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
int MyRank() const
Return the MPI rank within this object's communicator.
void Get2DPatchNets(const Vector &Nodes, int vdim)
void GetBdrElementTopo(Array< Element * > &boundary) const
void SetKnotsFromPatches()
void GetElementTopo(Array< Element * > &elements) const
const Array< int > & GetOrders() const
Read-only access to the orders of all knot vectors.
Array< int > f_spaceOffsets
Array< bool > activeBdrElem
int GetVDim() const
Returns vector dimension.
void LoadSolution(std::istream &input, GridFunction &sol) const
FiniteElementSpace * FESpace()
void SetElement(int e) const
int GetOrder() const
If all orders are identical, return that number. Otherwise, return NURBSFECollection::VariableOrder.
void Swap(Array< T > &, Array< T > &)
double p(const Vector &x, double t)
void LoadFE(int i, const FiniteElement *FE) const
const NURBSExtension * GetNURBSext() const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Array< int > v_meshOffsets
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void SetBdrPatchVertexMap(int p, const KnotVector *kv[], int *okv)
void GetBdrPatchKnotVectors(int p, Array< KnotVector * > &kv)
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Array< int > p_meshOffsets
KnotVector & operator=(const KnotVector &kv)
void SetPatchVertexMap(int p, const KnotVector *kv[])
Helper struct for defining a connectivity table, see Table::MakeFromList.
Table * GetElementDofTable()
void SetOrderFromOrders()
void CalcShape(Vector &shape, int i, double xi) const
void CalcDnShape(Vector &gradn, int n, int i, double xi) const
void PrintFunctions(std::ostream &out, int samples=11) const
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void GetPatchNets(const Vector &Nodes, int vdim)
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
NURBSPatch * Revolve3D(NURBSPatch &patch, double n[], double ang, int times)
Array< int > p_spaceOffsets
void SetSolutionVector(Vector &Nodes, int vdim)
void Get3DElementTopo(Array< Element * > &elements) const
void SetBdrPatchDofMap(int p, const KnotVector *kv[], int *okv)
void GenerateElementDofTable()
void Difference(const KnotVector &kv, Vector &diff) const
virtual void SetOrder() const
Update the NURBSFiniteElement according to the currently set knot vectors.
void Print(std::ostream &out) const
void DeleteAll()
Delete all dynamically allocated memory, resetting all dimensions to zero.
void CalcDShape(Vector &grad, int i, double xi) const
void DegreeElevate(int rel_degree, int degree=16)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void Get2DElementTopo(Array< Element * > &elements) const
Data type line segment element.
Array< int > bel_to_patch
void SetPatchDofMap(int p, const KnotVector *kv[])
void GetBdrElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of boundary element i.
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
void ConvertToPatches(const Vector &Nodes)
void KnotInsert(Array< KnotVector * > &kv)
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
void Get2DBdrElementTopo(Array< Element * > &boundary) const
int SetLoopDirection(int dir)