21void createGluingDataSpace(
const gsGeometry<real_t> & patch,
const gsBasis<real_t> & basis,
index_t dir,
22 gsBSplineBasis<real_t> & result,
index_t p_tilde,
index_t r_tilde);
24void createPlusSpace(
const gsGeometry<real_t> & patch, gsBasis<real_t> & basis,
index_t dir, gsBSplineBasis<real_t> & res_plus);
26void createMinusSpace(
const gsGeometry<real_t> & patch, gsBasis<real_t> & basis,
index_t dir, gsBSplineBasis<real_t> & res_minus);
28void createEdgeSpace(
const gsGeometry<real_t> & patch, gsBasis<real_t> & basis,
index_t dir, gsBSplineBasis<real_t> & basis_plus,
29 gsBSplineBasis<real_t> & basis_minus, gsBSplineBasis<real_t> & basis_gD, gsTensorBSplineBasis<2, real_t> & result);
31void createEdgeSpace(
const gsGeometry<real_t> & patch, gsBasis<real_t> & basis,
index_t dir, gsBSplineBasis<real_t> & basis_plus,
32 gsBSplineBasis<real_t> & basis_minus, gsTensorBSplineBasis<2, real_t> & result);
34void createVertexSpace(
const gsGeometry<real_t> & patch, gsBasis<real_t> & basis,
bool isInterface_1,
bool isInterface_2,
35 gsTensorBSplineBasis<2, real_t> & result,
index_t p_tilde,
index_t r_tilde);
44 mutable gsMapData<T> _tmp;
50 typedef memory::shared_ptr< gsAlpha > Ptr;
53 typedef memory::unique_ptr< gsAlpha > uPtr;
55 gsAlpha(gsGeometry<T> & geo,
index_t uv) :
56 _geo(geo), m_uv(uv), _alpha_piece(nullptr)
61 ~gsAlpha() {
delete _alpha_piece; }
63 GISMO_CLONE_FUNCTION(gsAlpha)
65 short_t domainDim()
const {
return 1;}
67 short_t targetDim()
const {
return 1;}
69 mutable gsAlpha<T> * _alpha_piece;
71 const gsFunction<T> & piece(
const index_t k)
const
75 _alpha_piece =
new gsAlpha(*
this);
80 void eval_into(
const gsMatrix<T>& u, gsMatrix<T>& result)
const
82 result.resize( targetDim() , u.cols() );
84 if(_geo.parDim() == _geo.targetDim())
86 gsMatrix<T> uv, ev, D0;
87 uv.setZero(2, u.cols());
92 for (
index_t i = 0; i < uv.cols(); i++) {
93 _geo.jacobian_into(uv.col(i), ev);
94 uv(0, i) = gamma * ev.determinant();
98 else if(_geo.parDim() + 1 == _geo.targetDim())
100 gsMatrix<T> uv, ev, D0, N, Duv;
101 uv.setZero(2, u.cols());
104 for (
index_t i = 0; i < uv.cols(); i++) {
106 _geo.deriv_into(uv.col(i), ev);
107 N.row(0) = ev.row(2) * ev.row(5) - ev.row(4) * ev.row(3);
108 N.row(1) = ev.row(4) * ev.row(1) - ev.row(0) * ev.row(5);
109 N.row(2) = ev.row(0) * ev.row(3) - ev.row(2) * ev.row(1);
115 Duv.row(0) = ev.row(0);
116 Duv.row(1) = ev.row(2);
117 Duv.row(2) = ev.row(4);
126 Duv.row(0) = ev.row(1);
127 Duv.row(1) = ev.row(3);
128 Duv.row(2) = ev.row(5);
135 uv(0, i) = 1.0 / Duv.norm() * gsEigen::numext::sqrt(D0.determinant());
150 gsGeometry<T> & _geo;
151 mutable gsMapData<T> _tmp;
157 typedef memory::shared_ptr< gsBeta > Ptr;
160 typedef memory::unique_ptr< gsBeta > uPtr;
162 gsBeta(gsGeometry<T> & geo,
index_t uv) :
163 _geo(geo), m_uv(uv), _beta_piece(nullptr)
168 ~gsBeta() {
delete _beta_piece; }
170 GISMO_CLONE_FUNCTION(gsBeta)
172 short_t domainDim()
const {
return 1;}
174 short_t targetDim()
const {
return 1;}
176 mutable gsBeta<T> * _beta_piece;
178 const gsFunction<T> & piece(
const index_t k)
const
182 _beta_piece =
new gsBeta(*
this);
187 void eval_into(
const gsMatrix<T>& u, gsMatrix<T>& result)
const
189 result.resize( targetDim() , u.cols() );
191 if(_geo.parDim() == _geo.targetDim())
193 gsMatrix<T> uv, ev, D0;
195 uv.setZero(2,u.cols());
200 for(
index_t i = 0; i < uv.cols(); i++)
202 _geo.jacobian_into(uv.col(i),ev);
204 T D1 = T(1.0)/ D0.norm();
205 uv(0,i) = - gamma * D1 * D1 * ev.col(1).transpose() * ev.col(0);
209 else if(_geo.parDim() + 1 == _geo.targetDim())
211 gsMatrix<T> uv, ev, D0;
213 uv.setZero(2,u.cols());
218 for(
index_t i = 0; i < uv.cols(); i++)
220 _geo.jacobian_into(uv.col(i),ev);
222 T D1 = T(1.0)/ D0.norm();
223 uv(0,i) = - gamma * D1 * D1 * ev.col(1).transpose() * ev.col(0);
238 gsGeometry<T> & _geo;
240 gsBSpline<T> m_basis_beta;
241 gsBSplineBasis<T> m_basis_plus;
243 gsBasis<T> & m_basis;
245 mutable gsMapData<T> _tmp;
253 typedef memory::shared_ptr< gsTraceBasis > Ptr;
256 typedef memory::unique_ptr< gsTraceBasis > uPtr;
258 gsTraceBasis(gsGeometry<T> & geo,
259 gsBSpline<T> basis_beta,
260 gsBSplineBasis<T> basis_plus,
265 _geo(geo), m_basis_beta(basis_beta), m_basis_plus(basis_plus), m_basis(
basis),
266 m_isboundary(isboundary), m_bfID(bfID), m_uv(uv), _traceBasis_piece(nullptr)
273 ~gsTraceBasis() {
delete _traceBasis_piece; }
275 GISMO_CLONE_FUNCTION(gsTraceBasis)
277 short_t domainDim()
const {
return 2;}
279 short_t targetDim()
const {
return 1;}
281 mutable gsTraceBasis<T> * _traceBasis_piece;
283 const gsFunction<T> & piece(
const index_t k)
const
287 _traceBasis_piece =
new gsTraceBasis(*
this);
288 return *_traceBasis_piece;
292 void eval_into(
const gsMatrix<T>& u, gsMatrix<T>& result)
const
294 result.resize( targetDim() , u.cols() );
297 gsBSplineBasis<T> bsp_temp =
dynamic_cast<gsBSplineBasis<T> &
>(m_basis.component(1-m_uv));
299 T p = bsp_temp.degree();
300 T tau_1 = bsp_temp.knots().at(p + 1);
302 gsMatrix<T> beta, N_0, N_1, N_i_plus, der_N_i_plus;
305 m_basis_beta.eval_into(u.row(m_uv),beta);
307 beta.setZero(1, u.cols());
309 m_basis.component(1-m_uv).evalSingle_into(0,u.row(1-m_uv),N_0);
310 m_basis.component(1-m_uv).evalSingle_into(1,u.row(1-m_uv),N_1);
312 m_basis_plus.evalSingle_into(m_bfID,u.row(m_uv),N_i_plus);
313 m_basis_plus.derivSingle_into(m_bfID,u.row(m_uv),der_N_i_plus);
315 gsMatrix<T> temp = beta.cwiseProduct(der_N_i_plus);
316 result = N_i_plus.cwiseProduct(N_0 + N_1) - temp.cwiseProduct(N_1) * tau_1 / p;
327 gsGeometry<T> & _geo;
329 gsBSpline<T> m_basis_alpha;
330 gsBSplineBasis<T> m_basis_minus;
332 gsBasis<T> & m_basis;
334 mutable gsMapData<T> _tmp;
342 typedef memory::shared_ptr< gsNormalDerivBasis > Ptr;
345 typedef memory::unique_ptr< gsNormalDerivBasis > uPtr;
347 gsNormalDerivBasis(gsGeometry<T> & geo,
348 gsBSpline<T> basis_alpha,
349 gsBSplineBasis<T> basis_minus,
354 _geo(geo), m_basis_alpha(basis_alpha), m_basis_minus(basis_minus), m_basis(
basis),
355 m_isboundary(isboundary), m_bfID(bfID), m_uv(uv), _normalDerivBasis_piece(nullptr)
361 ~gsNormalDerivBasis() {
delete _normalDerivBasis_piece; }
363 GISMO_CLONE_FUNCTION(gsNormalDerivBasis)
365 short_t domainDim()
const {
return 2;}
367 short_t targetDim()
const {
return 1;}
369 mutable gsNormalDerivBasis<T> * _normalDerivBasis_piece;
371 const gsFunction<T> & piece(
const index_t k)
const
375 _normalDerivBasis_piece =
new gsNormalDerivBasis(*
this);
376 return *_normalDerivBasis_piece;
380 void eval_into(
const gsMatrix<T>& u, gsMatrix<T>& result)
const
382 result.resize( targetDim() , u.cols() );
385 gsBSplineBasis<T> bsp_temp =
dynamic_cast<gsBSplineBasis<T> &
>(m_basis.component(1-m_uv));
387 T p = bsp_temp.degree();
388 T tau_1 = bsp_temp.knots().at(p + 1);
390 gsMatrix<T> alpha, N_1, N_j_minus;
393 m_basis_alpha.eval_into(u.row(m_uv),alpha);
395 alpha.setOnes(1, u.cols());
397 m_basis.component(1-m_uv).evalSingle_into(1,u.row(1-m_uv),N_1);
399 m_basis_minus.evalSingle_into(m_bfID,u.row(m_uv),N_j_minus);
402 result = (m_uv == 0 ? T(-1.0) : T(1.0)) * alpha.cwiseProduct(N_j_minus.cwiseProduct(N_1)) * tau_1 / p;
404 result = (m_uv == 0 ? T(-1.0) : T(1.0)) * alpha.cwiseProduct(N_j_minus.cwiseProduct(N_1));
415 const gsGeometry<T> & m_geo;
416 gsBasis<T> & m_basis;
418 std::vector<gsBSpline<T>> m_alpha;
419 std::vector<gsBSpline<T>> m_beta;
421 std::vector<gsBSplineBasis<T>> m_basis_plus;
422 std::vector<gsBSplineBasis<T>> m_basis_minus;
424 const gsMatrix<T> m_Phi;
425 const std::vector<bool> m_kindOfEdge;
429 mutable gsMapData<T> _tmp;
433 typedef memory::shared_ptr< gsVertexBasis > Ptr;
436 typedef memory::unique_ptr< gsVertexBasis > uPtr;
438 gsVertexBasis(
const gsGeometry<T> & geo,
440 std::vector<gsBSpline<T>> alpha,
441 std::vector<gsBSpline<T>> beta,
442 std::vector<gsBSplineBasis<T>> basis_plus,
443 std::vector<gsBSplineBasis<T>> basis_minus,
444 const gsMatrix<T> Phi,
445 const std::vector<bool> kindOfEdge,
447 ) : m_geo(geo), m_basis(
basis), m_alpha(alpha), m_beta(beta), m_basis_plus(basis_plus), m_basis_minus(basis_minus),
448 m_Phi(Phi), m_kindOfEdge(kindOfEdge), m_bfID(bfID),
449 _vertexBasis_piece(nullptr)
454 ~gsVertexBasis() {
delete _vertexBasis_piece; }
456 GISMO_CLONE_FUNCTION(gsVertexBasis)
458 short_t domainDim()
const {
return 2;}
460 short_t targetDim()
const {
return 1;}
462 mutable gsVertexBasis<T> * _vertexBasis_piece;
464 const gsFunction<T> & piece(
const index_t k)
const
468 _vertexBasis_piece =
new gsVertexBasis(*
this);
469 return *_vertexBasis_piece;
474 void eval_into(
const gsMatrix<T>& u, gsMatrix<T>& result)
const
476 result.resize( targetDim() , u.cols() );
484 std::vector<gsMatrix<T>> c_0, c_1;
485 std::vector<gsMatrix <T>> c_0_plus, c_1_plus, c_2_plus;
486 std::vector<gsMatrix <T>> c_0_plus_deriv, c_1_plus_deriv, c_2_plus_deriv;
487 std::vector<gsMatrix <T>> c_0_minus, c_1_minus;
488 for (
index_t i = 0; i < 2; i++)
490 gsMatrix<T> b_0, b_1;
491 gsMatrix<T> b_0_plus, b_1_plus, b_2_plus;
492 gsMatrix<T> b_0_plus_deriv, b_1_plus_deriv, b_2_plus_deriv;
493 gsMatrix<T> b_0_minus, b_1_minus;
501 m_basis.component(i).evalSingle_into(0, u.row(i),b_0);
502 m_basis.component(i).evalSingle_into(1, u.row(i),b_1);
504 m_basis_plus[i].evalSingle_into(0, u.row(i),b_0_plus);
505 m_basis_plus[i].evalSingle_into(1, u.row(i),b_1_plus);
506 m_basis_plus[i].evalSingle_into(2, u.row(i),b_2_plus);
508 m_basis_plus[i].derivSingle_into(0, u.row(i),b_0_plus_deriv);
509 m_basis_plus[i].derivSingle_into(1, u.row(i),b_1_plus_deriv);
510 m_basis_plus[i].derivSingle_into(2, u.row(i),b_2_plus_deriv);
512 m_basis_minus[i].evalSingle_into(0, u.row(i),b_0_minus);
513 m_basis_minus[i].evalSingle_into(1, u.row(i),b_1_minus);
515 gsMatrix<T> b_1_0, b_1_minus_0;
516 m_basis.component(i).derivSingle_into(1, zero.row(i),b_1_0);
517 m_basis_minus[i].derivSingle_into(1, zero.row(i),b_1_minus_0);
519 T factor_b_1 = 1.0/b_1_0(0,0);
520 c_0.push_back(b_0 + b_1);
521 c_1.push_back(factor_b_1 * b_1);
523 T factor_b_1_minus = 1.0/b_1_minus_0(0,0);
524 c_0_minus.push_back(b_0_minus + b_1_minus);
525 c_1_minus.push_back(factor_b_1_minus * b_1_minus);
527 gsMatrix<T> der_b_1_plus_0, der2_b_1_plus_0, der2_b_2_plus_0;
528 m_basis_plus[i].derivSingle_into(1, zero.row(i), der_b_1_plus_0);
529 m_basis_plus[i].deriv2Single_into(1, zero.row(i), der2_b_1_plus_0);
530 m_basis_plus[i].deriv2Single_into(2, zero.row(i), der2_b_2_plus_0);
532 T factor_c_1_plus = 1.0/der_b_1_plus_0(0,0);
533 T factor2_c_1_plus = -der2_b_1_plus_0(0,0)/(der_b_1_plus_0(0,0)*der2_b_2_plus_0(0,0));
534 T factor_c_2_plus = 1.0/der2_b_2_plus_0(0,0);
536 c_0_plus.push_back(b_0_plus + b_1_plus + b_2_plus);
537 c_1_plus.push_back(factor_c_1_plus * b_1_plus + factor2_c_1_plus * b_2_plus);
538 c_2_plus.push_back(factor_c_2_plus * b_2_plus );
540 c_0_plus_deriv.push_back(b_0_plus_deriv + b_1_plus_deriv + b_2_plus_deriv);
541 c_1_plus_deriv.push_back(factor_c_1_plus * b_1_plus_deriv + factor2_c_1_plus * b_2_plus_deriv);
542 c_2_plus_deriv.push_back(factor_c_2_plus * b_2_plus_deriv);
545 std::vector<gsMatrix<T>> alpha, beta, alpha_0, beta_0, alpha_deriv, beta_deriv;
546 gsMatrix < T > temp_mat;
549 m_alpha[0].eval_into(u.row(0),temp_mat);
550 alpha.push_back(temp_mat);
552 m_alpha[0].eval_into(zero.row(0),temp_mat);
553 alpha_0.push_back(temp_mat);
555 m_alpha[0].deriv_into(zero.row(0),temp_mat);
556 alpha_deriv.push_back(temp_mat);
558 m_beta[0].eval_into(u.row(0),temp_mat);
559 beta.push_back(temp_mat);
561 m_beta[0].eval_into(zero.row(0),temp_mat);
562 beta_0.push_back(temp_mat);
564 m_beta[0].deriv_into(zero.row(0),temp_mat);
565 beta_deriv.push_back(temp_mat);
569 temp_mat.setOnes(1, u.cols());
570 alpha.push_back(temp_mat);
572 temp_mat.setOnes(1, zero.cols());
573 alpha_0.push_back(temp_mat);
575 temp_mat.setZero(1, zero.cols());
576 alpha_deriv.push_back(temp_mat);
578 temp_mat.setZero(1, u.cols());
579 beta.push_back(temp_mat);
581 temp_mat.setZero(1, zero.cols());
582 beta_0.push_back(temp_mat);
584 temp_mat.setZero(1, zero.cols());
585 beta_deriv.push_back(temp_mat);
588 if (m_kindOfEdge[1]) {
589 m_alpha[1].eval_into(u.row(1), temp_mat);
590 alpha.push_back(temp_mat);
592 m_alpha[1].eval_into(zero.row(0), temp_mat);
593 alpha_0.push_back(temp_mat);
595 m_alpha[1].deriv_into(zero.row(0), temp_mat);
596 alpha_deriv.push_back(temp_mat);
598 m_beta[1].eval_into(u.row(1), temp_mat);
599 beta.push_back(temp_mat);
601 m_beta[1].eval_into(zero.row(0), temp_mat);
602 beta_0.push_back(temp_mat);
604 m_beta[1].deriv_into(zero.row(0), temp_mat);
605 beta_deriv.push_back(temp_mat);
609 temp_mat.setOnes(1, u.cols());
610 alpha.push_back(temp_mat);
612 temp_mat.setOnes(1, zero.cols());
613 alpha_0.push_back(temp_mat);
615 temp_mat.setZero(1, zero.cols());
616 alpha_deriv.push_back(temp_mat);
618 temp_mat.setZero(1, u.cols());
619 beta.push_back(temp_mat);
621 temp_mat.setZero(1, zero.cols());
622 beta_0.push_back(temp_mat);
624 temp_mat.setZero(1, zero.cols());
625 beta_deriv.push_back(temp_mat);
629 gsMatrix<T> geo_jac, geo_der2;
630 geo_jac = m_geo.jacobian(zero);
631 geo_der2 = m_geo.deriv2(zero);
636 gsMatrix<T> dd_ik_plus, dd_ik_minus;
637 gsMatrix<T> dd_ik_minus_deriv, dd_ik_plus_deriv;
638 dd_ik_minus = -1.0/(alpha_0[0](0,0)) * (geo_jac.col(1) +
639 beta_0[0](0,0) * geo_jac.col(0));
641 dd_ik_plus = 1.0/(alpha_0[1](0,0)) * (geo_jac.col(0) +
642 beta_0[1](0,0) * geo_jac.col(1));
644 gsMatrix<T> geo_deriv2_12(m_geo.targetDim(),1), geo_deriv2_11(m_geo.targetDim(),1), geo_deriv2_22(m_geo.targetDim(),1);
645 geo_deriv2_12.row(0) = geo_der2.row(2);
646 geo_deriv2_12.row(1) = geo_der2.row(5);
647 geo_deriv2_11.row(0) = geo_der2.row(0);
648 geo_deriv2_11.row(1) = geo_der2.row(3);
649 geo_deriv2_22.row(0) = geo_der2.row(1);
650 geo_deriv2_22.row(1) = geo_der2.row(4);
652 if(m_geo.parDim() + 1 == m_geo.targetDim())
654 geo_deriv2_12.row(2) = m_geo.deriv2(zero).row(8);
655 geo_deriv2_11.row(2) = m_geo.deriv2(zero).row(6);
656 geo_deriv2_22.row(2) = m_geo.deriv2(zero).row(7);
659 gsMatrix<T> alpha_squared_u = alpha_0[0]*alpha_0[0];
660 gsMatrix<T> alpha_squared_v = alpha_0[1]*alpha_0[1];
662 dd_ik_minus_deriv = -1.0/(alpha_squared_u(0,0)) *
663 ((geo_deriv2_12 + (beta_deriv[0](0,0) * geo_jac.col(0) +
664 beta_0[0](0,0) * geo_deriv2_11))*alpha_0[0](0,0) -
665 (geo_jac.col(1) + beta_0[0](0,0) * geo_jac.col(0)) *
666 alpha_deriv[0](0,0));
668 dd_ik_plus_deriv = 1.0/(alpha_squared_v(0,0)) *
669 ((geo_deriv2_12 + (beta_deriv[1](0,0) * geo_jac.col(1) +
670 beta_0[1](0,0) * geo_deriv2_22))*alpha_0[1](0,0) -
671 (geo_jac.col(0) + beta_0[1](0,0) * geo_jac.col(1)) *
672 alpha_deriv[1](0,0));
675 std::vector<gsMatrix<T>> d_ik;
676 d_ik.push_back(m_Phi.col(0));
677 d_ik.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * geo_jac.col(0) );
678 d_ik.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * geo_jac.col(1) );
680 if(m_geo.parDim() + 1 == m_geo.targetDim())
682 d_ik.push_back( (geo_jac(0,0) * m_Phi.row(4).transpose() + geo_jac(1,0) * m_Phi.row(7).transpose() + geo_jac(2,0) * m_Phi.row(10).transpose()) * geo_jac(0,1) +
683 (geo_jac(0,0) * m_Phi.row(5).transpose() + geo_jac(1,0) * m_Phi.row(8).transpose() + geo_jac(2,0) * m_Phi.row(11).transpose()) * geo_jac(1,1) +
684 (geo_jac(0,0) * m_Phi.row(6).transpose() + geo_jac(1,0) * m_Phi.row(9).transpose() + geo_jac(2,0) * m_Phi.row(12).transpose()) * geo_jac(2,1) +
685 m_Phi.block(1, 0, 1, 6).transpose() * geo_der2.row(2) +
686 m_Phi.block(2, 0, 1, 6).transpose() * geo_der2.row(5) +
687 m_Phi.block(3, 0, 1, 6).transpose() * geo_der2.row(8) );
692 d_ik.push_back((geo_jac(0, 0) * m_Phi.col(3) + geo_jac(1, 0) * m_Phi.col(4)) * geo_jac(0, 1) +
693 (geo_jac(0, 0) * m_Phi.col(4) + geo_jac(1, 0) * m_Phi.col(5)) * geo_jac(1, 1) +
694 m_Phi.block(0, 1, 6, 1) * geo_der2.row(2) +
695 m_Phi.block(0, 2, 6, 1) * geo_der2.row(5));
698 std::vector<gsMatrix<T>> d_ilik_minus, d_ilik_plus;
699 d_ilik_minus.push_back(m_Phi.col(0));
700 d_ilik_minus.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * geo_jac.col(0));
702 if(m_geo.parDim() + 1 == m_geo.targetDim())
704 d_ilik_minus.push_back( (geo_jac(0,0) * m_Phi.row(4).transpose() + geo_jac(1,0) * m_Phi.row(7).transpose() + geo_jac(2,0) * m_Phi.row(10).transpose()) * geo_jac(0,0) +
705 (geo_jac(0,0) * m_Phi.row(5).transpose() + geo_jac(1,0) * m_Phi.row(8).transpose() + geo_jac(2,0) * m_Phi.row(11).transpose()) * geo_jac(1,0) +
706 (geo_jac(0,0) * m_Phi.row(6).transpose() + geo_jac(1,0) * m_Phi.row(9).transpose() + geo_jac(2,0) * m_Phi.row(12).transpose()) * geo_jac(2,0) +
707 m_Phi.block(1, 0, 1, 6).transpose() * geo_der2.row(0) +
708 m_Phi.block(2, 0, 1, 6).transpose() * geo_der2.row(3) +
709 m_Phi.block(3, 0, 1, 6).transpose() * geo_der2.row(6) );
713 d_ilik_minus.push_back((geo_jac(0, 0) * m_Phi.col(3) + geo_jac(1, 0) * m_Phi.col(4)) * geo_jac(0, 0) +
714 (geo_jac(0, 0) * m_Phi.col(4) + geo_jac(1, 0) * m_Phi.col(5)) * geo_jac(1, 0) +
715 m_Phi.block(0, 1, 6, 1) * geo_der2.row(0) +
716 m_Phi.block(0, 2, 6, 1) * geo_der2.row(3));
719 d_ilik_minus.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * dd_ik_minus);
721 if(m_geo.parDim() + 1 == m_geo.targetDim())
723 d_ilik_minus.push_back( (geo_jac(0,0) * m_Phi.row(4).transpose() + geo_jac(1,0) * m_Phi.row(7).transpose() + geo_jac(2,0) * m_Phi.row(10).transpose()) * dd_ik_minus(0,0) +
724 (geo_jac(0,0) * m_Phi.row(5).transpose() + geo_jac(1,0) * m_Phi.row(8).transpose() + geo_jac(2,0) * m_Phi.row(11).transpose()) * dd_ik_minus(1,0) +
725 (geo_jac(0,0) * m_Phi.row(6).transpose() + geo_jac(1,0) * m_Phi.row(9).transpose() + geo_jac(2,0) * m_Phi.row(12).transpose()) * dd_ik_minus(2,0) +
726 m_Phi.block(1, 0, 1, 6).transpose() * dd_ik_minus_deriv.row(0) +
727 m_Phi.block(2, 0, 1, 6).transpose() * dd_ik_minus_deriv.row(1) +
728 m_Phi.block(3, 0, 1, 6).transpose() * dd_ik_minus_deriv.row(2) );
732 d_ilik_minus.push_back((geo_jac(0, 0) * m_Phi.col(3) + geo_jac(1, 0) * m_Phi.col(4)) * dd_ik_minus(0, 0) +
733 (geo_jac(0, 0) * m_Phi.col(4) + geo_jac(1, 0) * m_Phi.col(5)) * dd_ik_minus(1, 0) +
734 m_Phi.block(0, 1, 6, 1) * dd_ik_minus_deriv.row(0) +
735 m_Phi.block(0, 2, 6, 1) * dd_ik_minus_deriv.row(1));
738 d_ilik_plus.push_back(m_Phi.col(0));
739 d_ilik_plus.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * geo_jac.col(1));
741 if(m_geo.parDim() + 1 == m_geo.targetDim())
743 d_ilik_plus.push_back( (geo_jac(0,1) * m_Phi.row(4).transpose() + geo_jac(1,1) * m_Phi.row(7).transpose() + geo_jac(2,1) * m_Phi.row(10).transpose()) * geo_jac(0,1) +
744 (geo_jac(0,1) * m_Phi.row(5).transpose() + geo_jac(1,1) * m_Phi.row(8).transpose() + geo_jac(2,1) * m_Phi.row(11).transpose()) * geo_jac(1,1) +
745 (geo_jac(0,1) * m_Phi.row(6).transpose() + geo_jac(1,1) * m_Phi.row(9).transpose() + geo_jac(2,1) * m_Phi.row(12).transpose()) * geo_jac(2,1) +
746 m_Phi.block(1, 0, 1, 6).transpose() * geo_der2.row(1) +
747 m_Phi.block(2, 0, 1, 6).transpose() * geo_der2.row(4) +
748 m_Phi.block(3, 0, 1, 6).transpose() * geo_der2.row(7) );
752 d_ilik_plus.push_back((geo_jac(0, 1) * m_Phi.col(3) + geo_jac(1, 1) * m_Phi.col(4)) * geo_jac(0, 1) +
753 (geo_jac(0, 1) * m_Phi.col(4) + geo_jac(1, 1) * m_Phi.col(5)) * geo_jac(1, 1) +
754 m_Phi.block(0, 1, 6, 1) * geo_der2.row(1) +
755 m_Phi.block(0, 2, 6, 1) * geo_der2.row(4));
758 d_ilik_plus.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * dd_ik_plus);
760 if(m_geo.parDim() + 1 == m_geo.targetDim())
762 d_ilik_plus.push_back( (geo_jac(0,1) * m_Phi.row(4).transpose() + geo_jac(1,1) * m_Phi.row(7).transpose() + geo_jac(2,1) * m_Phi.row(10).transpose()) * dd_ik_plus(0,0) +
763 (geo_jac(0,1) * m_Phi.row(5).transpose() + geo_jac(1,1) * m_Phi.row(8).transpose() + geo_jac(2,1) * m_Phi.row(11).transpose()) * dd_ik_plus(1,0) +
764 (geo_jac(0,1) * m_Phi.row(6).transpose() + geo_jac(1,1) * m_Phi.row(9).transpose() + geo_jac(2,1) * m_Phi.row(12).transpose()) * dd_ik_plus(2,0) +
765 m_Phi.block(1, 0, 1, 6).transpose() * dd_ik_plus_deriv.row(0) +
766 m_Phi.block(2, 0, 1, 6).transpose() * dd_ik_plus_deriv.row(1) +
767 m_Phi.block(3, 0, 1, 6).transpose() * dd_ik_plus_deriv.row(2) );
771 d_ilik_plus.push_back((geo_jac(0, 1) * m_Phi.col(3) + geo_jac(1, 1) * m_Phi.col(4)) * dd_ik_plus(0, 0) +
772 (geo_jac(0, 1) * m_Phi.col(4) + geo_jac(1, 1) * m_Phi.col(5)) * dd_ik_plus(1, 0) +
773 m_Phi.block(0, 1, 6, 1) * dd_ik_plus_deriv.row(0) +
774 m_Phi.block(0, 2, 6, 1) * dd_ik_plus_deriv.row(1));
778 result = d_ilik_minus.at(0)(m_bfID,0) * (c_0_plus.at(0).cwiseProduct(c_0.at(1)) -
779 beta[0].cwiseProduct(c_0_plus_deriv.at(0).cwiseProduct(c_1.at(1)))) +
780 d_ilik_minus.at(1)(m_bfID,0) * (c_1_plus.at(0).cwiseProduct(c_0.at(1)) -
781 beta[0].cwiseProduct(c_1_plus_deriv.at(0).cwiseProduct(c_1.at(1)))) +
782 d_ilik_minus.at(2)(m_bfID,0) * (c_2_plus.at(0).cwiseProduct(c_0.at(1)) -
783 beta[0].cwiseProduct(c_2_plus_deriv.at(0).cwiseProduct(c_1.at(1)))) -
784 d_ilik_minus.at(3)(m_bfID,0) * alpha[0].cwiseProduct(c_0_minus.at(0).cwiseProduct(c_1.at(1))) -
785 d_ilik_minus.at(4)(m_bfID,0) * alpha[0].cwiseProduct(c_1_minus.at(0).cwiseProduct(c_1.at(1)));
791 result += d_ilik_plus.at(0)(m_bfID,0) * (c_0_plus.at(1).cwiseProduct(c_0.at(0)) -
792 beta[1].cwiseProduct(c_0_plus_deriv.at(1).cwiseProduct(c_1.at(0)))) +
793 d_ilik_plus.at(1)(m_bfID,0) * (c_1_plus.at(1).cwiseProduct(c_0.at(0)) -
794 beta[1].cwiseProduct(c_1_plus_deriv.at(1).cwiseProduct(c_1.at(0)))) +
795 d_ilik_plus.at(2)(m_bfID,0) * (c_2_plus.at(1).cwiseProduct(c_0.at(0)) -
796 beta[1].cwiseProduct(c_2_plus_deriv.at(1).cwiseProduct(c_1.at(0)))) +
797 d_ilik_plus.at(3)(m_bfID,0) * alpha[1].cwiseProduct(c_0_minus.at(1).cwiseProduct(c_1.at(0))) +
798 d_ilik_plus.at(4)(m_bfID,0) * alpha[1].cwiseProduct(c_1_minus.at(1).cwiseProduct(c_1.at(0)));
800 result -= d_ik.at(0)(m_bfID,0) * c_0.at(0).cwiseProduct(c_0.at(1)) + d_ik.at(2)(m_bfID,0) * c_0.at(0).cwiseProduct(c_1.at(1)) +
801 d_ik.at(1)(m_bfID,0) * c_1.at(0).cwiseProduct(c_0.at(1)) + d_ik.at(3)(m_bfID,0) * c_1.at(0).cwiseProduct(c_1.at(1));
const gsBasis< T > & basis(const index_t k) const
Helper which casts and returns the k-th piece of this function set as a gsBasis.
Definition gsFunctionSet.hpp:33
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
This is the interface of all objects that computes functions on points like gsBasis,...
The G+Smo namespace, containing all definitions for the library.
@ NEED_JACOBIAN
Jacobian of the object.
Definition gsForwardDeclarations.h:75