G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsApproxC1Utils.h
Go to the documentation of this file.
1
13#pragma once
14
15
17
18namespace gismo
19{
20
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);
23
24void createPlusSpace(const gsGeometry<real_t> & patch, gsBasis<real_t> & basis, index_t dir, gsBSplineBasis<real_t> & res_plus);
25
26void createMinusSpace(const gsGeometry<real_t> & patch, gsBasis<real_t> & basis, index_t dir, gsBSplineBasis<real_t> & res_minus);
27
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);
30
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);
33
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);
36
37// Input is parametric coordinates of 1-D \a mp
38template <class T>
39class gsAlpha : public gismo::gsFunction<T>
40{
41
42protected:
43 gsGeometry<T> & _geo;
44 mutable gsMapData<T> _tmp;
45 index_t m_uv;
46
47
48public:
50 typedef memory::shared_ptr< gsAlpha > Ptr;
51
53 typedef memory::unique_ptr< gsAlpha > uPtr;
54
55 gsAlpha(gsGeometry<T> & geo, index_t uv) :
56 _geo(geo), m_uv(uv), _alpha_piece(nullptr)
57 {
58 _tmp.flags = NEED_JACOBIAN;
59 }
60
61 ~gsAlpha() { delete _alpha_piece; }
62
63 GISMO_CLONE_FUNCTION(gsAlpha)
64
65 short_t domainDim() const {return 1;}
66
67 short_t targetDim() const {return 1;}
68
69 mutable gsAlpha<T> * _alpha_piece; // why do we need this?
70
71 const gsFunction<T> & piece(const index_t k) const
72 {
73 GISMO_UNUSED(k);
74 //delete _alpha_piece;
75 _alpha_piece = new gsAlpha(*this);
76 return *_alpha_piece;
77 }
78
79 // Input is parametric coordinates of 1-D \a mp
80 void eval_into(const gsMatrix<T>& u, gsMatrix<T>& result) const
81 {
82 result.resize( targetDim() , u.cols() );
83
84 if(_geo.parDim() == _geo.targetDim()) // Planar
85 {
86 gsMatrix<T> uv, ev, D0;
87 uv.setZero(2, u.cols());
88 uv.row(m_uv) = u; // u
89
90 T gamma = 1.0;
91
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();
95 }
96 result = uv.row(0);
97 }
98 else if(_geo.parDim() + 1 == _geo.targetDim()) // Surface
99 {
100 gsMatrix<T> uv, ev, D0, N, Duv;
101 uv.setZero(2, u.cols());
102 uv.row(m_uv) = u; // u
103
104 for (index_t i = 0; i < uv.cols(); i++) {
105 N.setZero(3,1);
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);
110
111 D0.setZero(3,3);
112 Duv.setZero(3,1);
113 if (m_uv == 0) // u
114 {
115 Duv.row(0) = ev.row(0);
116 Duv.row(1) = ev.row(2);
117 Duv.row(2) = ev.row(4);
118 D0.col(0) = Duv;
119 D0.col(2) = N;
120 D0(0,1) = ev(1,0);
121 D0(1,1) = ev(3,0);
122 D0(2,1) = ev(5,0);
123 }
124 else if (m_uv == 1) // v
125 {
126 Duv.row(0) = ev.row(1);
127 Duv.row(1) = ev.row(3);
128 Duv.row(2) = ev.row(5);
129 D0.col(1) = Duv;
130 D0.col(2) = N;
131 D0(0,0) = ev(0,0);
132 D0(1,0) = ev(2,0);
133 D0(2,0) = ev(4,0);
134 }
135 uv(0, i) = 1.0 / Duv.norm() * gsEigen::numext::sqrt(D0.determinant());
136 }
137 result = uv.row(0);
138 //gsDebugVar(result);
139 //gsDebugVar(m_uv);
140 }
141 }
142};
143
144
145template <class T>
146class gsBeta : public gismo::gsFunction<T>
147{
148
149protected:
150 gsGeometry<T> & _geo;
151 mutable gsMapData<T> _tmp;
152 index_t m_uv;
153
154
155public:
157 typedef memory::shared_ptr< gsBeta > Ptr;
158
160 typedef memory::unique_ptr< gsBeta > uPtr;
161
162 gsBeta(gsGeometry<T> & geo, index_t uv) :
163 _geo(geo), m_uv(uv), _beta_piece(nullptr)
164 {
165 _tmp.flags = NEED_JACOBIAN;
166 }
167
168 ~gsBeta() { delete _beta_piece; }
169
170 GISMO_CLONE_FUNCTION(gsBeta)
171
172 short_t domainDim() const {return 1;}
173
174 short_t targetDim() const {return 1;}
175
176 mutable gsBeta<T> * _beta_piece; // why do we need this?
177
178 const gsFunction<T> & piece(const index_t k) const
179 {
180 GISMO_UNUSED(k);
181 //delete _beta_piece;
182 _beta_piece = new gsBeta(*this);
183 return *_beta_piece;
184 }
185
186 // Input is parametric coordinates of 1-D \a mp
187 void eval_into(const gsMatrix<T>& u, gsMatrix<T>& result) const
188 {
189 result.resize( targetDim() , u.cols() );
190
191 if(_geo.parDim() == _geo.targetDim()) // Planar
192 {
193 gsMatrix<T> uv, ev, D0;
194
195 uv.setZero(2,u.cols());
196 uv.row(m_uv) = u; // u
197
198 T gamma = 1.0;
199
200 for(index_t i = 0; i < uv.cols(); i++)
201 {
202 _geo.jacobian_into(uv.col(i),ev);
203 D0 = ev.col(m_uv);
204 T D1 = T(1.0)/ D0.norm();
205 uv(0,i) = - gamma * D1 * D1 * ev.col(1).transpose() * ev.col(0);
206 }
207 result = uv.row(0);
208 }
209 else if(_geo.parDim() + 1 == _geo.targetDim()) // Surface
210 {
211 gsMatrix<T> uv, ev, D0;
212
213 uv.setZero(2,u.cols());
214 uv.row(m_uv) = u; // u
215
216 T gamma = 1.0;
217
218 for(index_t i = 0; i < uv.cols(); i++)
219 {
220 _geo.jacobian_into(uv.col(i),ev);
221 D0 = ev.col(m_uv);
222 T D1 = T(1.0)/ D0.norm();
223 uv(0,i) = - gamma * D1 * D1 * ev.col(1).transpose() * ev.col(0);
224 }
225 result = uv.row(0);
226 //gsDebugVar(result);
227 //gsDebugVar(m_uv);
228 }
229 }
230
231};
232
233template <class T>
234class gsTraceBasis : public gismo::gsFunction<T>
235{
236
237protected:
238 gsGeometry<T> & _geo;
239
240 gsBSpline<T> m_basis_beta;
241 gsBSplineBasis<T> m_basis_plus;
242
243 gsBasis<T> & m_basis;
244
245 mutable gsMapData<T> _tmp;
246
247 bool m_isboundary;
248 const index_t m_bfID, m_uv;
249
250
251public:
253 typedef memory::shared_ptr< gsTraceBasis > Ptr;
254
256 typedef memory::unique_ptr< gsTraceBasis > uPtr;
257
258 gsTraceBasis(gsGeometry<T> & geo,
259 gsBSpline<T> basis_beta,
260 gsBSplineBasis<T> basis_plus,
261 gsBasis<T> & basis,
262 bool isboundary,
263 const index_t bfID,
264 const index_t uv) :
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)
267 {
268 //_tmp.flags = NEED_JACOBIAN;
269
270 //createPlusSpace(geo, basis, m_uv, m_basis_plus); // Not efficient
271 }
272
273 ~gsTraceBasis() { delete _traceBasis_piece; }
274
275 GISMO_CLONE_FUNCTION(gsTraceBasis)
276
277 short_t domainDim() const {return 2;}
278
279 short_t targetDim() const {return 1;}
280
281 mutable gsTraceBasis<T> * _traceBasis_piece; // why do we need this?
282
283 const gsFunction<T> & piece(const index_t k) const
284 {
285 GISMO_UNUSED(k);
286 //delete _traceBasis_piece;
287 _traceBasis_piece = new gsTraceBasis(*this);
288 return *_traceBasis_piece;
289 }
290
291 // Input is parametric coordinates of 2-D \a mp
292 void eval_into(const gsMatrix<T>& u, gsMatrix<T>& result) const
293 {
294 result.resize( targetDim() , u.cols() );
295
296 // tau/p
297 gsBSplineBasis<T> bsp_temp = dynamic_cast<gsBSplineBasis<T> & >(m_basis.component(1-m_uv));
298
299 T p = bsp_temp.degree();
300 T tau_1 = bsp_temp.knots().at(p + 1); // p + 2
301
302 gsMatrix<T> beta, N_0, N_1, N_i_plus, der_N_i_plus;
303
304 if (!m_isboundary)
305 m_basis_beta.eval_into(u.row(m_uv),beta); // 1-dir == PatchID
306 else
307 beta.setZero(1, u.cols());
308
309 m_basis.component(1-m_uv).evalSingle_into(0,u.row(1-m_uv),N_0); // u
310 m_basis.component(1-m_uv).evalSingle_into(1,u.row(1-m_uv),N_1); // u
311
312 m_basis_plus.evalSingle_into(m_bfID,u.row(m_uv),N_i_plus); // v
313 m_basis_plus.derivSingle_into(m_bfID,u.row(m_uv),der_N_i_plus);
314
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;
317 }
318
319};
320
321
322template <class T>
323class gsNormalDerivBasis : public gismo::gsFunction<T>
324{
325
326protected:
327 gsGeometry<T> & _geo;
328
329 gsBSpline<T> m_basis_alpha;
330 gsBSplineBasis<T> m_basis_minus;
331
332 gsBasis<T> & m_basis;
333
334 mutable gsMapData<T> _tmp;
335
336 bool m_isboundary;
337 const index_t m_bfID, m_uv;
338
339
340public:
342 typedef memory::shared_ptr< gsNormalDerivBasis > Ptr;
343
345 typedef memory::unique_ptr< gsNormalDerivBasis > uPtr;
346
347 gsNormalDerivBasis(gsGeometry<T> & geo,
348 gsBSpline<T> basis_alpha,
349 gsBSplineBasis<T> basis_minus,
350 gsBasis<T> & basis,
351 bool isboundary,
352 const index_t bfID,
353 const index_t uv) :
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)
356 {
357 //_tmp.flags = NEED_JACOBIAN;
358 //createMinusSpace(geo, basis, m_uv, m_basis_minus);
359 }
360
361 ~gsNormalDerivBasis() { delete _normalDerivBasis_piece; }
362
363 GISMO_CLONE_FUNCTION(gsNormalDerivBasis)
364
365 short_t domainDim() const {return 2;}
366
367 short_t targetDim() const {return 1;}
368
369 mutable gsNormalDerivBasis<T> * _normalDerivBasis_piece; // why do we need this?
370
371 const gsFunction<T> & piece(const index_t k) const
372 {
373 GISMO_UNUSED(k);
374 //delete _normalDerivBasis_piece;
375 _normalDerivBasis_piece = new gsNormalDerivBasis(*this);
376 return *_normalDerivBasis_piece;
377 }
378
379 // Input is parametric coordinates of 2-D \a mp
380 void eval_into(const gsMatrix<T>& u, gsMatrix<T>& result) const
381 {
382 result.resize( targetDim() , u.cols() );
383
384 // tau/p
385 gsBSplineBasis<T> bsp_temp = dynamic_cast<gsBSplineBasis<T> & >(m_basis.component(1-m_uv));
386
387 T p = bsp_temp.degree();
388 T tau_1 = bsp_temp.knots().at(p + 1); // p + 2
389
390 gsMatrix<T> alpha, N_1, N_j_minus;
391
392 if (!m_isboundary)
393 m_basis_alpha.eval_into(u.row(m_uv),alpha); // 1-dir == PatchID
394 else
395 alpha.setOnes(1, u.cols());
396
397 m_basis.component(1-m_uv).evalSingle_into(1,u.row(1-m_uv),N_1); // u
398
399 m_basis_minus.evalSingle_into(m_bfID,u.row(m_uv),N_j_minus); // v
400
401 if (!m_isboundary)
402 result = (m_uv == 0 ? T(-1.0) : T(1.0)) * alpha.cwiseProduct(N_j_minus.cwiseProduct(N_1)) * tau_1 / p;
403 else
404 result = (m_uv == 0 ? T(-1.0) : T(1.0)) * alpha.cwiseProduct(N_j_minus.cwiseProduct(N_1));
405 }
406
407};
408
409
410template <class T>
411class gsVertexBasis : public gismo::gsFunction<T>
412{
413
414protected:
415 const gsGeometry<T> & m_geo;
416 gsBasis<T> & m_basis;
417
418 std::vector<gsBSpline<T>> m_alpha;
419 std::vector<gsBSpline<T>> m_beta;
420
421 std::vector<gsBSplineBasis<T>> m_basis_plus;
422 std::vector<gsBSplineBasis<T>> m_basis_minus;
423
424 const gsMatrix<T> m_Phi;
425 const std::vector<bool> m_kindOfEdge;
426
427 const index_t m_bfID;
428
429 mutable gsMapData<T> _tmp;
430
431public:
433 typedef memory::shared_ptr< gsVertexBasis > Ptr;
434
436 typedef memory::unique_ptr< gsVertexBasis > uPtr;
437
438 gsVertexBasis(const gsGeometry<T> & geo,
439 gsBasis<T> & basis,
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,
446 const index_t bfID
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)
450 {
451
452 }
453
454 ~gsVertexBasis() { delete _vertexBasis_piece; }
455
456 GISMO_CLONE_FUNCTION(gsVertexBasis)
457
458 short_t domainDim() const {return 2;}
459
460 short_t targetDim() const {return 1;}
461
462 mutable gsVertexBasis<T> * _vertexBasis_piece; // why do we need this?
463
464 const gsFunction<T> & piece(const index_t k) const
465 {
466 GISMO_UNUSED(k);
467 //delete _vertexBasis_piece;
468 _vertexBasis_piece = new gsVertexBasis(*this);
469 return *_vertexBasis_piece;
470 }
471
472 // Input is parametric coordinates of 2-D \a mp
473 // TODO IMPROVE THE FUNCTION
474 void eval_into(const gsMatrix<T>& u, gsMatrix<T>& result) const
475 {
476 result.resize( targetDim() , u.cols() );
477 result.setZero();
478
479 // Computing c, c+ and c-
480 // Point zero
481 gsMatrix<T> zero;
482 zero.setZero(2,1);
483
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++) // i == 0 == u , i == 1 == v
489 {
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;
494
495 //gsBSplineBasis<T> bsp_temp = dynamic_cast<gsBSplineBasis<T> & >(m_basis.component(i));
496 //gsBSplineBasis<T> bsp_temp_2 = dynamic_cast<gsBSplineBasis<T> & >(m_basis.component(i));
497 //T p = bsp_temp.degree();
498 //T h_geo = bsp_temp.knots().at(p + 1);
499 //T h_geo_2 = bsp_temp_2.knots().at(p + 1);
500
501 m_basis.component(i).evalSingle_into(0, u.row(i),b_0); // first
502 m_basis.component(i).evalSingle_into(1, u.row(i),b_1); // second
503
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);
507
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);
511
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);
514
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);
518
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);
522
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);
526
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);
531
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);
535
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 );
539
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);
543 }
544
545 std::vector<gsMatrix<T>> alpha, beta, alpha_0, beta_0, alpha_deriv, beta_deriv;
546 gsMatrix < T > temp_mat;
547 if (m_kindOfEdge[0])
548 {
549 m_alpha[0].eval_into(u.row(0),temp_mat); // 1-dir == PatchID
550 alpha.push_back(temp_mat); // u
551
552 m_alpha[0].eval_into(zero.row(0),temp_mat); // 1-dir == PatchID
553 alpha_0.push_back(temp_mat); // u
554
555 m_alpha[0].deriv_into(zero.row(0),temp_mat); // 1-dir == PatchID
556 alpha_deriv.push_back(temp_mat); // u
557
558 m_beta[0].eval_into(u.row(0),temp_mat); // 1-dir == PatchID
559 beta.push_back(temp_mat); // u
560
561 m_beta[0].eval_into(zero.row(0),temp_mat); // 1-dir == PatchID
562 beta_0.push_back(temp_mat); // u
563
564 m_beta[0].deriv_into(zero.row(0),temp_mat); // 1-dir == PatchID
565 beta_deriv.push_back(temp_mat); // u
566 }
567 else
568 {
569 temp_mat.setOnes(1, u.cols());
570 alpha.push_back(temp_mat); // u
571
572 temp_mat.setOnes(1, zero.cols());
573 alpha_0.push_back(temp_mat); // u
574
575 temp_mat.setZero(1, zero.cols());
576 alpha_deriv.push_back(temp_mat); // u
577
578 temp_mat.setZero(1, u.cols());
579 beta.push_back(temp_mat); // u
580
581 temp_mat.setZero(1, zero.cols());
582 beta_0.push_back(temp_mat); // u
583
584 temp_mat.setZero(1, zero.cols());
585 beta_deriv.push_back(temp_mat); // u
586 }
587
588 if (m_kindOfEdge[1]) {
589 m_alpha[1].eval_into(u.row(1), temp_mat); // 1-dir == PatchID
590 alpha.push_back(temp_mat); // v
591
592 m_alpha[1].eval_into(zero.row(0), temp_mat); // 1-dir == PatchID
593 alpha_0.push_back(temp_mat); // v
594
595 m_alpha[1].deriv_into(zero.row(0), temp_mat); // 1-dir == PatchID
596 alpha_deriv.push_back(temp_mat); // v
597
598 m_beta[1].eval_into(u.row(1), temp_mat); // 1-dir == PatchID
599 beta.push_back(temp_mat); // v
600
601 m_beta[1].eval_into(zero.row(0), temp_mat); // 1-dir == PatchID
602 beta_0.push_back(temp_mat); // v
603
604 m_beta[1].deriv_into(zero.row(0), temp_mat); // 1-dir == PatchID
605 beta_deriv.push_back(temp_mat); // v
606 }
607 else
608 {
609 temp_mat.setOnes(1, u.cols());
610 alpha.push_back(temp_mat); // u
611
612 temp_mat.setOnes(1, zero.cols());
613 alpha_0.push_back(temp_mat); // u
614
615 temp_mat.setZero(1, zero.cols());
616 alpha_deriv.push_back(temp_mat); // u
617
618 temp_mat.setZero(1, u.cols());
619 beta.push_back(temp_mat); // u
620
621 temp_mat.setZero(1, zero.cols());
622 beta_0.push_back(temp_mat); // u
623
624 temp_mat.setZero(1, zero.cols());
625 beta_deriv.push_back(temp_mat); // u
626 }
627
628 // Geo data:
629 gsMatrix<T> geo_jac, geo_der2;
630 geo_jac = m_geo.jacobian(zero);
631 geo_der2 = m_geo.deriv2(zero);
632
633
634
635 // Compute dd^^(i_k) and dd^^(i_k-1)
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));
640
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));
643
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);
651
652 if(m_geo.parDim() + 1 == m_geo.targetDim()) // Surface
653 {
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);
657 }
658
659 gsMatrix<T> alpha_squared_u = alpha_0[0]*alpha_0[0];
660 gsMatrix<T> alpha_squared_v = alpha_0[1]*alpha_0[1];
661
662 dd_ik_minus_deriv = -1.0/(alpha_squared_u(0,0)) * // N^2
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));
667
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));
673
674 // Comupute d_(0,0)^(i_k), d_(1,0)^(i_k), d_(0,1)^(i_k), d_(1,1)^(i_k) ; i_k == 2
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) ); // deriv into u
678 d_ik.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * geo_jac.col(1) ); // deriv into v
679
680 if(m_geo.parDim() + 1 == m_geo.targetDim()) // In the surface case the dimension of the second derivative vector is 9x1
681 {
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) );
688
689 }
690 else
691 {
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)); // Hessian
696 }
697 // Compute d_(*,*)^(il,ik)
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));
701
702 if(m_geo.parDim() + 1 == m_geo.targetDim()) // In the surface case the dimension of the second derivative vector is 9x1
703 {
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) );
710 }
711 else
712 {
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));
717 }
718
719 d_ilik_minus.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * dd_ik_minus);
720
721 if(m_geo.parDim() + 1 == m_geo.targetDim()) // In the surface case the dimension of the second derivative vector is 9x1
722 {
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) );
729 }
730 else
731 {
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));
736 }
737
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));
740
741 if(m_geo.parDim() + 1 == m_geo.targetDim()) // In the surface case the dimension of the second derivative vector is 9x1
742 {
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) );
749 }
750 else
751 {
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));
756 }
757
758 d_ilik_plus.push_back(m_Phi.block(1, 0, m_geo.targetDim(), 6).transpose() * dd_ik_plus);
759
760 if(m_geo.parDim() + 1 == m_geo.targetDim()) // In the surface case the dimension of the second derivative vector is 9x1
761 {
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) );
768 }
769 else
770 {
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));
775 }
776
777
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))); // f*_(ik-1,ik)
786
787 //if (kindOfEdge[0])
788 //rhsVals.at(i).setZero();
789
790 //if (!kindOfEdge[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))); // f*_(ik+1,ik)
799
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)); // f*_(ik)
802
803
804 }
805
806};
807
808}
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