G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsBarrierCore.hpp
Go to the documentation of this file.
1
19#pragma once
20
22
23namespace gismo {
24
25template<short_t d, typename T>
27 gsOptionList options;
28
29 // Print output 0: no print, 1: summary, 2: iterations and summary
30 options.addInt("Verbose", "Print output level", 0);
31
32 // Initialization Method: 0 Coons' patch (default), 1 Spring patch, 2: Cross-Ap. patch
33 options.addInt("InitialMethod", "Initialization method", 0);
34
35 // Parameterization Method: 0 Barrier patch (default), 1 Penalty patch, 2: PDE patch
36 options.addInt("ParamMethod", "Parameterization method", 0);
37
38 // Parameter and stopping criteria for foldover elimination step
39 options.addReal("ff_Delta", "Delta for foldover-free optimization", 0.05);
40 options.addInt("ff_MaxIterations",
41 "Max iterations for quality improvement",
42 1e4);
43 options.addReal("ff_MinGradientLength",
44 "Min gradient length for foldover-free optimization",
45 1e-20);
46 options.addReal("ff_MinStepLength",
47 "Min step length for foldover-free optimization",
48 1e-20);
49
50 // Parameters and stopping criteria for quality improvement step
51 options.addInt("qi_MaxIterations",
52 "Max iterations for quality improvement",
53 1e4);
54 options.addReal("qi_MinGradientLength",
55 "Min gradient length for quality improvement",
56 1e-4);
57 options.addReal("qi_MinStepLength",
58 "Min step length for quality improvement",
59 1e-4);
60
61 // Set quadrature rules for objective function and gradient evaluation
62 options.addReal("quA",
63 "Quadrature points: quA*deg + quB; For patchRule: Order of target space",
64 1.0);
65 options.addInt("quB",
66 "Quadrature points: quA*deg + quB; For patchRule: Regularity of target space",
67 1);
68 options.addInt("quRule",
69 "Quadrature rule [1:GaussLegendre,2:GaussLobatto,3:PatchRule]",
70 1);
71 options.addInt("overInt", "Apply over-integration?", 0);
72
73 // Preconditioner type for AA: 0: NO preconditioning, 1: Full Jacobian preconditioner, 2: Diagonal Jacobian preconditioner, 3: Diagonal Block Jacobian preconditioner
74 options.addInt("AAPreconditionType", "Preconditioner type for AA", 0);
75
76 // Update the preconditioner every N_update steps
77 options.addInt("N_update", "update preconditioner every N_update steps", 10);
78
79 // window size m for our preconditioned AA solver
80 options.addInt("AAwindowsize", "window size for preconditioned AA solver", 5);
81
82 // need improve the parameterization quality for PDE-based method?
83 options.addSwitch("needPDEH1", "improve quality by H1 discrezation?", true);
84
85 return options;
86}
87
89template<short_t d, typename T>
91 return computeAreaInterior(mp);
92}
93
95template<short_t d, typename T>
97 // Creating a multi-basis from the provided multi-patch
98 gsMultiBasis<T> multiBasis(multiPatch);
99
100 // Initializing an expression evaluator
101 gsExprEvaluator<T> evaluator;
102
103 // Setting integration elements
104 evaluator.setIntegrationElements(multiBasis);
105
106 // Getting the geometry map
107 geometryMap geometry = evaluator.getMap(multiPatch);
108
109 // Computing the area by integrating over the determinant of the Jacobian of the geometry map
110 T area = evaluator.integral(jac(geometry).det());
111
112 return area;
113}
114
115template<short_t d, typename T>
117{
118 GISMO_ERROR("NO_IMPLEMENTATION");
119}
120
121enum class ParamMethod {
122 PDEPatch,
123 PenaltyPatch,
124 PenaltyPatch2,
125 BarrierPatch,
126 VariationalHarmonicPatch,
127 // Add new methods here...
128};
129
130template<short_t d, typename T>
132 const gsDofMapper &mapper,
133 const gsOptionList &options) {
134 gsMultiPatch<T> result;
135
136 ParamMethod
137 method = static_cast<ParamMethod>(options.askInt("ParamMethod", 0));
138
139 switch (method) {
140 case ParamMethod::PenaltyPatch: {
141#ifdef gsHLBFGS_ENABLED
142 result = computePenaltyPatch(mp, mapper, options);
143#else
144 GISMO_ERROR("PenaltyPatch not available without gsHLBFGS. Please "
145 "compile with gsHLBGFS enabled");
146#endif
147 break;
148 }
149 case ParamMethod::PenaltyPatch2: {
150#ifdef gsHLBFGS_ENABLED
151 result = computePenaltyPatch2(mp, mapper, options);
152#else
153 GISMO_ERROR("PenaltyPatch2 not available without gsHLBFGS. Please "
154 "compile with gsHLBGFS enabled");
155#endif
156 break;
157 }
158 case ParamMethod::BarrierPatch: {
159#ifdef gsHLBFGS_ENABLED
160 result = computeBarrierPatch(mp, mapper, options);
161#else
162 GISMO_ERROR("BarrierPatch not available without gsHLBFGS. Please "
163 "compile with gsHLBGFS enabled");
164#endif
165 break;
166 }
167 case ParamMethod::VariationalHarmonicPatch: {
168#ifdef gsHLBFGS_ENABLED
169 result = computeVHPatch(mp, mapper, options);
170#else
171 GISMO_ERROR("VHPatch not available without gsHLBFGS. Please "
172 "compile with gsHLBGFS enabled");
173#endif
174 break;
175 }
176 case ParamMethod::PDEPatch:
177 default: {
178 if (method != ParamMethod::PDEPatch) {
179 gsWarn << "Invalid ParamMethod value. Defaulting to PDEPatch.\n";
180 }
181 result = computePDEPatch(mp, mapper, options);
182 break;
183 }
184 }
185 return result;
186}
187
188#ifdef gsHLBFGS_ENABLED
189// Modified Variational Harmonic Method
190template<short_t d, typename T>
193 const gismo::gsDofMapper &mapper,
194 const gismo::gsOptionList &options) {
195
196 verboseLog("Start variational harmonic parameterization construction...",
197 options.askInt("Verbose", 0));
198
199 // get initial guess vector
200 gsVector<T> initialGuessVector = convertMultiPatchToFreeVector(mp, mapper);
201
202 gsObjVHPt<d, T> objVHPt(mp, mapper);
203 objVHPt.applyOptions(options);
204
205 gsHLBFGS<T> optimizer(&objVHPt);
206 setOptimizerOptions<T>(optimizer, options);
207
208 optimizer.solve(initialGuessVector);
209 gsMultiPatch<T> result = mp;
210 convertFreeVectorToMultiPatch<T>(optimizer.currentDesign(), mapper, result);
211
212 verboseLog("Finished!", options.askInt("Verbose", 0));
213
214 return result;
215}
216
217template<short_t d, typename T>
218gsMultiPatch<T>
219gsBarrierCore<d, T>::computePenaltyPatch(const gsMultiPatch<T> &mp,
220 const gsDofMapper &mapper,
221 const gsOptionList &options) {
222
223 verboseLog("Start penalty function-based parameterization construction...",
224 options.askInt("Verbose", 0));
225
226 // Compute scaledArea and initial guess vector
227 T scaledArea = computeAreaInterior(mp);
228
229 // get initial guess vector
230 gsVector<T> initialGuessVector = convertMultiPatchToFreeVector(mp, mapper);
231
232 gsObjPenaltyPt<d, T> objPenaltyPt(mp, mapper);
233 gsOptionList thisOptions = options;
234 thisOptions.addReal("qi_lambda1", "Sets the lambda_1 value for Emips", 1.0);
235 thisOptions.addReal("qi_lambda2",
236 "Sets the lambda 2 value for Eunif",
237 1.0 / pow(scaledArea, 2));
238 objPenaltyPt.applyOptions(thisOptions);
239
240 gsHLBFGS<T> optimizer(&objPenaltyPt);
241 setOptimizerOptions<T>(optimizer, options);
242
243 optimizer.solve(initialGuessVector);
244 gsMultiPatch<T> result = mp;
245 convertFreeVectorToMultiPatch<T>(optimizer.currentDesign(), mapper, result);
246
247 verboseLog("Finished!", options.askInt("Verbose", 0));
248
249 return result;
250}
251
252template<short_t d, typename T>
253gsMultiPatch<T>
254gsBarrierCore<d, T>::computePenaltyPatch2(const gsMultiPatch<T> &mp,
255 const gsDofMapper &mapper,
256 const gsOptionList &options) {
257
258 verboseLog("Penalty function-based (2) parameterization construction ...",
259 options.askInt("Verbose", 0));
260
261 // Compute scaledArea and initial guess vector
262 T scaledArea = computeAreaInterior(mp);
263
264 // get initial guess vector
265 gsVector<T> initialGuessVector = convertMultiPatchToFreeVector(mp, mapper);
266
267 gsObjPenaltyPt2<d, T> objPenaltyPt(mp, mapper);
268 objPenaltyPt.setEps(options.askReal("penaltyFactor", 0.01) * scaledArea);
269
270 gsOptionList thisOptions = options;
271 thisOptions.addReal("qi_lambda1", "Sets the lambda_1 value for Emips", 1.0);
272 thisOptions.addReal("qi_lambda2",
273 "Sets the lambda 2 value for Eunif",
274 1.0 / pow(scaledArea, 2));
275 objPenaltyPt.applyOptions(thisOptions);
276
277 gsHLBFGS<T> optimizer(&objPenaltyPt);
278 setOptimizerOptions<T>(optimizer, options);
279
280 optimizer.solve(initialGuessVector);
281 gsMultiPatch<T> result = mp;
282 convertFreeVectorToMultiPatch<T>(optimizer.currentDesign(), mapper, result);
283
284 verboseLog("Finished!", options.askInt("Verbose", 0));
285
286 return result;
287}
288
289template<short_t d, typename T>
290gsMultiPatch<T>
291 gsBarrierCore<d,T>::computeBarrierPatch(const gsMultiPatch<T> &mp,
292 const gsDofMapper &mapper,
293 const gsOptionList &options) {
294
295 // Compute scaledArea and initial guess vector
296 T scaledArea = computeAreaInterior(mp);
297
298 // get initial guess vector
299 gsVector<T> initialGuessVector = convertMultiPatchToFreeVector(mp, mapper);
300
301 // STEP 2: foldover elimination step
302 foldoverElimination(mp, mapper, initialGuessVector, scaledArea, options);
303
304 // STEP 3: parameterization quality improvement
305 qualityImprovement(mp, mapper, initialGuessVector, scaledArea, options);
306
307 // Update the result with optimized parameterization
308 gsMultiPatch<T> result = mp;
309 convertFreeVectorToMultiPatch(initialGuessVector, mapper, result);
310
311 return result;
312}
313
314template<short_t d, typename T>
315void gsBarrierCore<d, T>::foldoverElimination(const gsMultiPatch<T> &mp,
316 const gsDofMapper &mapper,
317 gsVector<T> &initialGuessVector,
318 const T &scaledArea,
319 const gsOptionList &options) {
320
321 verboseLog("Start foldover elimination step...",
322 options.askInt("Verbose", 0));
323 constexpr T EPSILON = 1e-20;
324 constexpr int MAX_ITER = 10;
325 gsObjFoldoverFree<d, T> objFoldoverFree(mp, mapper);
326 objFoldoverFree.addOptions(options);
327
328 gsHLBFGS<T> optFoldoverFree(&objFoldoverFree);
329 optFoldoverFree.options().setInt("MaxIterations",
330 options.askInt("ff_MaxIterations", 1e4));
331 optFoldoverFree.options().setReal("MinGradientLength",
332 options.askReal("ff_MinGradientLength",
333 1e-12));
334 optFoldoverFree.options().setReal("MinStepLength",
335 options.askReal("ff_MinStepLength", 1e-12));
336 optFoldoverFree.options().setInt("Verbose", options.askInt("Verbose", 0));
337
338 T Efoldover = std::numeric_limits<T>::max();
339 for (index_t it = 0; it < MAX_ITER; ++it) {
340 T delta = pow(0.1, it) * 5e-2 * scaledArea; // parameter delta
341 objFoldoverFree.setDelta(delta);
342
343 optFoldoverFree.solve(initialGuessVector);
344
345 Efoldover = optFoldoverFree.objective();
346 initialGuessVector = optFoldoverFree.currentDesign();
347
348 if (Efoldover <= EPSILON) { break; }
349 }
350
351 if (Efoldover > EPSILON) {
352 throw std::runtime_error(
353 "Maximum iterations reached. The foldover-energy value is " +
354 std::to_string(Efoldover) +
355 ". This suggests there may be issues with the input data.");
356 }
357}
358
359template<short_t d, typename T>
360void gsBarrierCore<d, T>::qualityImprovement(const gsMultiPatch<T> &mp,
361 const gsDofMapper &mapper,
362 gsVector<T> &initialGuessVector,
363 const T &scaledArea,
364 const gsOptionList &options) {
365 verboseLog("Start parameterization quality improvement step...",
366 options.askInt("Verbose", 0));
367 gsObjQualityImprovePt<d, T> objQualityImprovePt(mp, mapper);
368
369 gsOptionList thisOptions = options;
370 thisOptions.addReal("qi_lambda1", "Sets the lambda_1 value for Emips", 1.0);
371 thisOptions.addReal("qi_lambda2",
372 "Sets the lambda 2 value for Eunif",
373 1.0 / pow(scaledArea, 2));
374 objQualityImprovePt.applyOptions(thisOptions);
375
376 gsHLBFGS<T> optQualityImprovePt(&objQualityImprovePt);
377 setOptimizerOptions<T>(optQualityImprovePt, options);
378
379 optQualityImprovePt.solve(initialGuessVector);
380 initialGuessVector = optQualityImprovePt.currentDesign();
381}
382
383template<short_t d, typename T>
384gsObjFoldoverFree<d, T>::gsObjFoldoverFree(const gsMultiPatch<T> &patches,
385 gsDofMapper mapper)
386 :
387 m_mp(patches),
388 m_mapper(std::move(mapper)),
389 m_mb(m_mp) {
391 m_assembler.setIntegrationElements(m_mb);
392 m_evaluator = gsExprEvaluator<T>(m_assembler);
393}
394
395template<short_t d, typename T>
396void gsObjFoldoverFree<d, T>::defaultOptions() {
397 // @Ye, make this reasonable default options
398 m_options.addReal("ff_Delta", "Sets the delta value", 1e-2);
399}
400
401template<short_t d, typename T>
402void gsObjFoldoverFree<d, T>::addOptions(const gsOptionList &options) {
403 m_options.update(options, gsOptionList::addIfUnknown);
404}
405
406template<short_t d, typename T>
407void gsObjFoldoverFree<d, T>::applyOptions(const gsOptionList &options) {
408 m_eps = m_options.getReal("ff_Delta");
409 m_evaluator.options().update(m_options, gsOptionList::addIfUnknown);
410}
411
412template<short_t d, typename T>
413T gsObjFoldoverFree<d, T>::evalObj(const gsAsConstVector<T> &u) const {
414 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
415 geometryMap G = m_evaluator.getMap(m_mp);
416
417 auto EfoldoverFree = (m_eps - jac(G).det()).ppartval();
418 return m_evaluator.integral(EfoldoverFree);
419}
420
421template<short_t d, typename T>
422void gsObjFoldoverFree<d, T>::gradObj_into(const gsAsConstVector<T> &u,
423 gsAsVector<T> &result) const {
424 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
425 geometryMap G = m_assembler.getMap(m_mp);
426
427 // Only call these once if their results don't change.
428 space space1 = m_assembler.getSpace(m_mb, d);
429 space1.setupMapper(m_mapper);
430
431 // |J|' w.r.t. physical coordinates x and y
432 auto derJacDet = frprod2(space1, jac(G).tr().adj());
433
434 gsConstantFunction<T> zeroFunc(gsVector<T>::Zero(d), d);
435 auto zeroVar = m_evaluator.getVariable(zeroFunc);
436
437 auto Eder = ternary(m_eps - jac(G).det(), -derJacDet, space1 * zeroVar);
438
439 m_assembler.initSystem();
440 m_assembler.assemble(Eder);
441
442 result.resize(m_assembler.rhs().rows());
443 std::copy(m_assembler.rhs().data(),
444 m_assembler.rhs().data() + m_assembler.rhs().rows(),
445 result.data());
446}
447
448template<short_t d, typename T>
449gsObjQualityImprovePt<d, T>::gsObjQualityImprovePt(
450 const gsMultiPatch<T> &patches,
451 gsDofMapper mapper)
452 :
453 m_mp(patches),
454 m_mapper(std::move(mapper)),
455 m_mb(m_mp) {
456 m_assembler.setIntegrationElements(m_mb);
457 m_evaluator = gsExprEvaluator<T>(m_assembler);
458// defaultOptions();
459}
460
461template<short_t d, typename T>
462void gsObjQualityImprovePt<d, T>::defaultOptions() {
463 // @Ye, make this reasonable default options
464 m_options.addReal("qi_lambda1", "Sets the lambda 1 value", 1.0);
465 m_options.addReal("qi_lambda2", "Sets the lambda 2 value", 1.0);
466}
467
468template<short_t d, typename T>
469void gsObjQualityImprovePt<d, T>::addOptions(const gsOptionList &options) {
470 m_options.update(options, gsOptionList::addIfUnknown);
471}
472
473template<short_t d, typename T>
474void gsObjQualityImprovePt<d, T>::applyOptions(const gsOptionList &options) {
475 m_options.update(options, gsOptionList::addIfUnknown);
476 m_lambda1 = m_options.getReal("qi_lambda1");
477 m_lambda2 = m_options.getReal("qi_lambda2");
478 m_evaluator.options().update(m_options, gsOptionList::addIfUnknown);
479}
480
481template<short_t d, typename T>
482T gsObjQualityImprovePt<d, T>::evalObj(const gsAsConstVector<T> &u) const {
483 return evalObj_impl<d>(u);
484}
485
486template<short_t d, typename T>
487template<short_t _d>
488typename std::enable_if<_d == 2, T>::type
489gsObjQualityImprovePt<d, T>::evalObj_impl(const gsAsConstVector<T> &u) const {
490 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
491 geometryMap G = m_evaluator.getMap(m_mp);
492
493 if (m_evaluator.min(jac(G).det()) < 0) {
494 return std::numeric_limits<T>::max();
495 } else {
496 auto Ewinslow = jac(G).sqNorm() / jac(G).det();
497 auto Euniform = pow(jac(G).det(), 2);
498
499 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
500 }
501}
502
503template<short_t d, typename T>
504template<short_t _d>
505typename std::enable_if<_d == 3, T>::type
506gsObjQualityImprovePt<d, T>::evalObj_impl(const gsAsConstVector<T> &u) const {
507// gsMultiPatch<T> mp = m_mp;
508 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
509
510 geometryMap G = m_evaluator.getMap(m_mp);
511
512 if (m_evaluator.min(jac(G).det()) < 0) {
513 return std::numeric_limits<T>::max();
514 } else {
515 auto Euniform = pow(jac(G).det(), 2);
516 auto Ewinslow = 0.5 * (jac(G).sqNorm() * jac(G).inv().sqNorm());
517
518// // another objective function term - Jiao et al. 2011
519// auto Ewinslow = jac(G).sqNorm() / pow(jac(G).det(), 2.0 / 3.0);
520
521 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
522 }
523}
524
525template<short_t d, typename T>
526void gsObjQualityImprovePt<d, T>::gradObj_into(const gsAsConstVector<T> &u,
527 gsAsVector<T> &result) const {
528 gradObj_into_impl<d>(u, result);
529}
530
531template<short_t d, typename T>
532template<short_t _d>
533typename std::enable_if<_d == 2, T>::type
534gsObjQualityImprovePt<d, T>::gradObj_into_impl(const gsAsConstVector<T> &u,
535 gsAsVector<T> &result) const {
536 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
537
538 geometryMap G = m_assembler.getMap(m_mp);
539
540 space space1 = m_assembler.getSpace(m_mb, d); // 1D space!!
541 space1.setupMapper(m_mapper);
542
543 // |J|' w.r.t. physical coordinates x and y
544 auto derJacDet = frprod2(space1, jac(G).tr().adj());
545
546 auto Ewinslow = jac(G).sqNorm() / jac(G).det();
547 auto derEwinslow = 2.0 / jac(G).det() * (frprod2(space1, jac(G))) -
548 Ewinslow.val() / jac(G).det() * derJacDet;
549 auto derEuniform = 2.0 * jac(G).det() * derJacDet;
550
551 m_assembler.initSystem();
552 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
553 result = gsAsVector<T>(const_cast<T *>(m_assembler.rhs().data()),
554 m_assembler.rhs().rows());
555 return EXIT_SUCCESS;
556}
557
558template<short_t d, typename T>
559template<short_t _d>
560typename std::enable_if<_d == 3, T>::type
561gsObjQualityImprovePt<d, T>::gradObj_into_impl(const gsAsConstVector<T> &u,
562 gsAsVector<T> &result) const {
563 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
564 geometryMap G = m_assembler.getMap(m_mp);
565
566 space space1 = m_assembler.getSpace(m_mb, d); // 1D space!!
567 space1.setupMapper(m_mapper);
568
569 // |J|' w.r.t. physical coordinates x and y
570 auto derJacDet = frprod2(space1, jac(G).tr().adj());
571
572 auto derEwinslow = frprod2(space1, (jac(G).inv().sqNorm() * jac(G) - jac(G)
573 .sqNorm() * (jac(G).tr() * jac(G) * jac(G).tr()).inv()));
574
575// // gradient of another objective function term
576// auto Ewinslow = jac(G).sqNorm() / pow(jac(G).det(), 2.0 / 3.0);
577// auto derEwinslow =
578// 2.0 * frprod2(space1, jac(G)) / pow(jac(G).det(), 2.0 / 3.0) -
579// 2.0 * Ewinslow.val() * derJacDet / jac(G).det() / 3.0;
580 auto derEuniform = 2.0 * jac(G).det() * derJacDet;
581
582 m_assembler.initSystem();
583 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
584 result = gsAsVector<T>(const_cast<T *>(m_assembler.rhs().data()),
585 m_assembler.rhs().rows());
586 return EXIT_SUCCESS;
587}
588
589template<short_t d, typename T>
590gsObjVHPt<d, T>::gsObjVHPt(const gsMultiPatch<T> &patches,
591 gsDofMapper mapper)
592 :
593 m_mp(patches),
594 m_mapper(std::move(mapper)),
595 m_mb(m_mp) {
597 m_assembler.setIntegrationElements(m_mb);
598 m_evaluator = gsExprEvaluator<T>(m_assembler);
599}
600
601template<short_t d, typename T>
602void gsObjVHPt<d, T>::defaultOptions() {
603 // @Ye, make this reasonable default options
604 m_options.addReal("qi_lambda1", "Sets the lambda 1 value", 1.0);
605 m_options.addReal("qi_lambda2", "Sets the lambda 2 value", 1.0);
606}
607
608template<short_t d, typename T>
609void gsObjVHPt<d, T>::addOptions(const gsOptionList &options) {
610 m_options.update(options, gsOptionList::addIfUnknown);
611}
612
613template<short_t d, typename T>
614void gsObjVHPt<d, T>::applyOptions(const gsOptionList &options) {
615 m_options.update(options, gsOptionList::addIfUnknown);
616 m_evaluator.options().update(m_options, gsOptionList::addIfUnknown);
617}
618
619template<short_t d, typename T>
620T gsObjVHPt<d, T>::evalObj(const gsAsConstVector<T> &u) const {
621 return evalObj_impl<d>(u);
622}
623
624template<short_t d, typename T>
625template<short_t _d>
626typename std::enable_if<_d == 2, T>::type
627gsObjVHPt<d, T>::evalObj_impl(const gsAsConstVector<T> &u) const {
628 // Convert the free vector to multi-patch
629 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
630
631 // Get the map for the geometry
632 geometryMap G = m_evaluator.getMap(m_mp);
633
634 // Get and set up the 1D space
635 space space1 = m_assembler.getSpace(m_mb, d);
636 space1.setupMapper(m_mapper);
637
638 // Compute the metric matrix
639 auto metricMat = jac(G).tr() * jac(G);
640
641 // Compute the Hessian
642 auto hessG = hess(G);
643
644 // Calculate the scale factor
645 auto scale = metricMat.trace();
646
647 // Calculate the Lx value
648 auto Lx = hessG % metricMat.adj() / scale.val();
649
650 // Build the nonlinear system
651 auto nonlinearSystem = frprod3(space1, Lx);
652
653 // Initialize and assemble the system
654 m_assembler.initSystem();
655 m_assembler.assemble(nonlinearSystem);
656
657 // Create the result vector
658 gsVector<T> result = gsAsVector<T>(
659 const_cast<T *>(m_assembler.rhs().data()),
660 m_assembler.rhs().rows());
661
662 return result.squaredNorm();
663}
664
665template<short_t d, typename T>
666template<short_t _d>
667typename std::enable_if<_d == 3, T>::type
668gsObjVHPt<d, T>::evalObj_impl(const gsAsConstVector<T> &u) const {
670}
671
672template<short_t d, typename T>
673void gsObjVHPt<d, T>::gradObj2_into(const gsAsConstVector<T> &u,
674 gsAsVector<T> &result) const {
675 gradObj_into_impl<d>(u, result);
676}
677
678template<short_t d, typename T>
679gsObjPenaltyPt<d, T>::gsObjPenaltyPt(const gsMultiPatch<T> &patches,
680 gsDofMapper mapper)
681 :
682 m_mp(patches),
683 m_mapper(std::move(mapper)),
684 m_mb(m_mp) {
686 m_assembler.setIntegrationElements(m_mb);
687 m_evaluator = gsExprEvaluator<T>(m_assembler);
688}
689
690template<short_t d, typename T>
691void gsObjPenaltyPt<d, T>::defaultOptions() {
692 // @Ye, make this reasonable default options
693 m_options.addReal("qi_lambda1", "Sets the lambda 1 value", 1.0);
694 m_options.addReal("qi_lambda2", "Sets the lambda 2 value", 1.0);
695}
696
697template<short_t d, typename T>
698void gsObjPenaltyPt<d, T>::addOptions(const gsOptionList &options) {
699 m_options.update(options, gsOptionList::addIfUnknown);
700}
701
702template<short_t d, typename T>
703void gsObjPenaltyPt<d, T>::applyOptions(const gsOptionList &options) {
704 m_options.update(options, gsOptionList::addIfUnknown);
705 m_lambda1 = m_options.getReal("qi_lambda1");
706 m_lambda2 = m_options.getReal("qi_lambda2");
707 m_evaluator.options().update(m_options, gsOptionList::addIfUnknown);
708}
709
710template<short_t d, typename T>
711T gsObjPenaltyPt<d, T>::evalObj(const gsAsConstVector<T> &u) const {
712 // Convert the free vector to multipatch
713 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
714
715 // Generate map for the geometry
716 geometryMap G = m_evaluator.getMap(m_mp);
717
718 // Generate epsilon value
719 gsConstantFunction<T> eps1(pow(m_eps, 2.0), d);
720 auto eps = m_evaluator.getVariable(eps1);
721
722 // Calculate chi value - penalty function
723 auto
724 chi = 0.5 * (jac(G).det() + pow(eps.val() + pow(jac(G).det(), 2.0), 0.5));
725
726 // Calculate Ewinslow and Euniform
727 auto Ewinslow = jac(G).sqNorm() / pow(chi, (2.0 / static_cast<T>(d)));
728 auto Euniform = pow(jac(G).det(), 2.0);
729
730 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
731// return evalObj_impl<d>(u);
732}
733
734template<short_t d, typename T>
735void gsObjPenaltyPt<d, T>::gradObj_into(const gsAsConstVector<T> &u,
736 gsAsVector<T> &result) const {
737 gradObj_into_impl<d>(u, result);
738}
739
740template<short_t d, typename T>
741template<short_t _d>
742typename std::enable_if<_d == 2, T>::type
743gsObjPenaltyPt<d, T>::gradObj_into_impl(const gsAsConstVector<T> &u,
744 gsAsVector<T> &result) const {
745
746 // Convert the free vector to multipatch
747 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
748
749 // Get map and space
750 geometryMap G = m_assembler.getMap(m_mp);
751 space space1 = m_assembler.getSpace(m_mb, d);
752 space1.setupMapper(m_mapper);
753
754 // Derivative of Jacobian determinant
755 auto derJacDet = jac(space1) % jac(G).tr().adj();
756
757 // Define constants for evaluation
758 gsConstantFunction<T> eps1(pow(m_eps, 2.0), d);
759 gsConstantFunction<T> unit1(1.0, d);
760 auto eps = m_evaluator.getVariable(eps1);
761 auto unit = m_evaluator.getVariable(unit1);
762
763 // Compute common term, chi and chip
764 auto commonTerm = pow(eps.val() + pow(jac(G).det(), 2.0), 0.5);
765 auto chi = 0.5 * (jac(G).det() + commonTerm);
766 auto chip = 0.5 * (unit.val() + jac(G).det() / commonTerm);
767
768 // Compute Ewinslow
769 auto Ewinslow = jac(G).sqNorm() / chi;
770
771 // Compute derivatives of Ewinslow and Euniform
772 auto derEwinslow = 1.0 / chi * (2.0 * frprod2(space1, jac(G)) -
773 Ewinslow.val() * chip * derJacDet);
774 auto derEuniform = 2.0 * jac(G).det() * derJacDet;
775
776 // Assemble the system
777 m_assembler.initSystem();
778 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
779
780 // Update the result
781 result = gsAsVector<T>(const_cast<T *>(m_assembler.rhs().data()),
782 m_assembler.rhs().rows());
783
784 return EXIT_SUCCESS;
785}
786
787template<short_t d, typename T>
788template<short_t _d>
789typename std::enable_if<_d == 3, T>::type
790gsObjPenaltyPt<d, T>::gradObj_into_impl(const gsAsConstVector<T> &u,
791 gsAsVector<T> &result) const {
792 const T twoThirds = 2.0 / 3.0;
793
794 // Convert the free vector to multipatch
795 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
796
797 // Get map and space
798 geometryMap G = m_assembler.getMap(m_mp);
799 space space1 = m_assembler.getSpace(m_mb, d); // 1D space!!
800 space1.setupMapper(m_mapper);
801
802 // Derivative of Jacobian determinant
803 auto derJacDet = frprod2(space1, jac(G).tr().adj());
804
805 // Define constants for evaluation
806 gsConstantFunction<T> eps1(pow(m_eps, 2.0), d);
807 gsConstantFunction<T> unit1(1.0, d);
808 auto eps = m_evaluator.getVariable(eps1);
809 auto unit = m_evaluator.getVariable(unit1);
810
811 // Compute common term, chi and chip
812 auto commonTerm = pow(eps.val() + pow(jac(G).det(), 2.0), 0.5);
813 auto chi = 0.5 * (jac(G).det() + commonTerm);
814 auto chip = 0.5 * (unit.val() + jac(G).det() / commonTerm);
815
816 // Compute Ewinslow
817 auto Ewinslow = jac(G).sqNorm() / pow(chi, twoThirds);
818
819 // Compute derivatives of Ewinslow and Euniform
820 auto derEwinslow = 2.0 * frprod2(space1, jac(G)) / pow(chi, twoThirds) -
821 twoThirds * Ewinslow / chi * chip * derJacDet;
822 auto derEuniform = 2.0 * jac(G).det() * derJacDet;
823
824 // Assemble the system
825 m_assembler.initSystem();
826 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
827
828 // Update the result
829 result = gsAsVector<T>(const_cast<T *>(m_assembler.rhs().data()),
830 m_assembler.rhs().rows());
831
832 return EXIT_SUCCESS;
833}
834
835// gsObjPenaltyPt2: my penalty function
836template<short_t d, typename T>
837gsObjPenaltyPt2<d, T>::gsObjPenaltyPt2(const gsMultiPatch<T> &patches,
838 gsDofMapper mapper)
839 :
840 m_mp(patches),
841 m_mapper(std::move(mapper)),
842 m_mb(m_mp) {
844 m_assembler.setIntegrationElements(m_mb);
845 m_evaluator = gsExprEvaluator<T>(m_assembler);
846}
847
848template<short_t d, typename T>
849void gsObjPenaltyPt2<d, T>::defaultOptions() {
850 // @Ye, make this reasonable default options
851 m_options.addReal("qi_lambda1", "Sets the lambda 1 value", 1.0);
852 m_options.addReal("qi_lambda2", "Sets the lambda 2 value", 1.0);
853}
854
855template<short_t d, typename T>
856void gsObjPenaltyPt2<d, T>::addOptions(const gsOptionList &options) {
857 m_options.update(options, gsOptionList::addIfUnknown);
858}
859
860template<short_t d, typename T>
861void gsObjPenaltyPt2<d, T>::applyOptions(const gsOptionList &options) {
862 m_options.update(options, gsOptionList::addIfUnknown);
863 m_lambda1 = m_options.getReal("qi_lambda1");
864 m_lambda2 = m_options.getReal("qi_lambda2");
865 m_evaluator.options().update(m_options, gsOptionList::addIfUnknown);
866}
867
868template<short_t d, typename T>
869T gsObjPenaltyPt2<d, T>::evalObj(const gsAsConstVector<T> &u) const {
870 return evalObj_impl<d>(u);
871}
872
873template<short_t d, typename T>
874template<short_t _d>
875typename std::enable_if<_d == 2, T>::type
876gsObjPenaltyPt2<d, T>::evalObj_impl(const gsAsConstVector<T> &u) const {
877 // Convert the free vector to multipatch
878 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
879
880 // Get the map G
881 geometryMap G = m_evaluator.getMap(m_mp);
882
883 // Constant function for evaluation
884 gsConstantFunction<T> epsilon(m_eps, d);
885 auto eps = m_evaluator.getVariable(epsilon);
886
887 // Calculation of chi with ternary operation
888 auto chiPPart = eps * ((jac(G).det() - eps.val()).exp());
889 auto chi =
890 ternary(eps.val() - jac(G).det(), chiPPart.val(), jac(G).det().val());
891
892 // Calculation of Ewinslow and Euniform
893 auto Ewinslow = jac(G).sqNorm() / chi;
894// auto Euniform = chi + 1.0 / chi;
895 auto Euniform = pow(jac(G).det(), 2.0);
896
897 // Evaluate the integral and return
898 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
899}
900
901template<short_t d, typename T>
902template<short_t _d>
903typename std::enable_if<_d == 3, T>::type
904gsObjPenaltyPt2<d, T>::evalObj_impl(const gsAsConstVector<T> &u) const {
905 // Convert the free vector to multipatch
906 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
907
908 // Get the map G
909 geometryMap G = m_evaluator.getMap(m_mp);
910
911 // Define epsilon constant function
912 gsConstantFunction<T> epsilon(m_eps, d);
913 auto eps = m_evaluator.getVariable(epsilon);
914
915 // Compute chi part and chi
916 auto chiPPart = eps * ((jac(G).det() - eps.val()).exp());
917 auto chi =
918 ternary(eps.val() - jac(G).det(), chiPPart.val(), jac(G).det().val());
919
920 // Compute Ewinslow
921 auto Ewinslow = 0.5 * (jac(G).sqNorm() * jac(G).inv().sqNorm()) *
922 pow(jac(G).det(), 2.0) / pow(chi, 2.0);
923
924 // Compute Euniform
925// auto Euniform = chi + 1.0 / chi;
926 auto Euniform = pow(jac(G).det(), 2.0);
927
928 // Compute and return the integral
929 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
930}
931
932template<short_t d, typename T>
933void gsObjPenaltyPt2<d, T>::gradObj_into(const gsAsConstVector<T> &u,
934 gsAsVector<T> &result) const {
935 gradObj_into_impl<d>(u, result);
936}
937
938template<short_t d, typename T>
939template<short_t _d>
940typename std::enable_if<_d == 2, T>::type
941gsObjPenaltyPt2<d, T>::gradObj_into_impl(const gsAsConstVector<T> &u,
942 gsAsVector<T> &result) const {
943 // Convert the free vector to multipatch
944 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
945
946 // Get the map G
947 geometryMap G = m_assembler.getMap(m_mp);
948
949 // Set up the space mapper
950 space spaceMapper = m_assembler.getSpace(m_mb, d);
951 spaceMapper.setupMapper(m_mapper);
952
953 // Compute the derivative of |J| with respect to physical coordinates x and y
954 auto derJacDet = frprod2(spaceMapper, jac(G).tr().adj());
955
956 // Define constant functions for epsilon and unity
957 gsConstantFunction<T> epsilon(m_eps, d);
958 gsConstantFunction<T> unity(1.0, d);
959
960 // Get the variables for epsilon and unity
961 auto eps = m_evaluator.getVariable(epsilon);
962 auto unit = m_evaluator.getVariable(unity);
963
964 // Define chi and its derivative chip
965 auto chiPPart = eps * ((jac(G).det() - eps.val()).exp());
966 auto chi =
967 ternary(eps.val() - jac(G).det(), chiPPart.val(), jac(G).det().val());
968 auto chip = ternary(eps.val() - jac(G).det(), chiPPart.val(), unit.val());
969
970 // Define Ewinslow and its derivative
971 auto Ewinslow = jac(G).sqNorm() / chi;
972 auto derEwinslow = 1.0 / chi * (2.0 * frprod2(spaceMapper, jac(G)) -
973 Ewinslow.val() * chip * derJacDet);
974
975 // Define the derivative of Euniform
976// auto derEuniform = (chip - chip / pow(chi, 2)) * derJacDet;
977 auto derEuniform = 2.0 * jac(G).det() * derJacDet;
978
979 // Assemble the system and set the result
980 m_assembler.initSystem();
981 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
982 result = gsAsVector<T>(const_cast<T *>(m_assembler.rhs().data()),
983 m_assembler.rhs().rows());
984
985 // Return success
986 return EXIT_SUCCESS;
987}
988
989template<short_t d, typename T>
990template<short_t _d>
991typename std::enable_if<_d == 3, T>::type
992gsObjPenaltyPt2<d, T>::gradObj_into_impl(const gsAsConstVector<T> &u,
993 gsAsVector<T> &result) const {
994 // Convert the free vector to multipatch
995 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
996
997 // Initialize the map G
998 geometryMap G = m_assembler.getMap(m_mp);
999
1000 // Initialize and setup the space
1001 space space1 = m_assembler.getSpace(m_mb, d);
1002 space1.setupMapper(m_mapper);
1003
1004 // Compute the derivative of |J| w.r.t. physical coordinates x and y
1005 auto derJacDet = frprod2(space1, jac(G).tr().adj());
1006
1007 // Initialize constant functions
1008 gsConstantFunction<T> epsilonFunction(m_eps, d);
1009 gsConstantFunction<T> unityFunction(1.0, d);
1010
1011 // Get the corresponding variables
1012 auto eps = m_evaluator.getVariable(epsilonFunction);
1013 auto unit = m_evaluator.getVariable(unityFunction);
1014
1015 // Compute the chi part
1016 auto chiPPart = eps * ((jac(G).det() - eps.val()).exp());
1017
1018 // Ternary operation to compute chi and chip
1019 auto chi =
1020 ternary(eps.val() - jac(G).det(), chiPPart.val(), jac(G).det().val());
1021 auto chip = ternary(eps.val() - jac(G).det(), chiPPart.val(), unit.val());
1022
1023 // Define additional computations
1024 auto jacFrobNorm2 = jac(G).sqNorm() * jac(G).inv().sqNorm();
1025 auto derJacFrob2 = frprod2(space1,
1026 (jac(G).inv().sqNorm() * jac(G) - jac(G).sqNorm()
1027 * (jac(G).tr() * jac(G) * jac(G).tr()).inv()));
1028
1029 auto derEwinslow = derJacFrob2 * pow(jac(G).det(), 2) / pow(chi, 2)
1030 + jacFrobNorm2 / pow(chi, 2)
1031 * (jac(G).det() - chip / chi * pow(jac(G).det(), 2)) * derJacDet;
1032// auto derEuniform = (chip - chip / pow(chi, 2)) * derJacDet;
1033 auto derEuniform = 2.0 * jac(G).det() * derJacDet;
1034
1035 // Initialize and assemble the system
1036 m_assembler.initSystem();
1037 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
1038
1039 // Compute and return the result
1040 result = gsAsVector<T>(const_cast<T *>(m_assembler.rhs().data()),
1041 m_assembler.rhs().rows());
1042 return EXIT_SUCCESS;
1043}
1044#endif
1045
1046template<short_t d, typename T>
1047gsMultiPatch<T>
1049 const gismo::gsDofMapper &mapper,
1050 const gismo::gsOptionList &options) {
1051
1052 // get initial guess vector
1053 gsVector<T> initialGuessVector = convertMultiPatchToFreeVector(mp, mapper);
1054
1055// std::string solverSetting = "Parameters setting:\n";
1056// solverSetting += "\t\t\t Window size = " +
1057// std::to_string(options.getInt("AAwindowsize")) + "\n";
1058// solverSetting += "\t\t\t Update preconditioner every " +
1059// std::to_string(options.getInt("N_update")) + " steps \n";
1060// solverSetting += "\t\t\t Preconditioner type: " +
1061// std::to_string(options.getInt("AAPreconditionType")) + "\n";
1062// solverSetting += "\t\t\t Need improve by H1? " +
1063// std::to_string(options.getSwitch("needPDEH1")) + "\n";
1064// verboseLog(solverSetting, options.askInt("Verbose", 0));
1065
1066 verboseLog("PDE-based parameterization construction ...\n",
1067 options.askInt("Verbose", 0));
1068
1069 gsExprEvaluator<T> evaluator;
1070 gsExprAssembler<T> assembler;
1071
1072 gsMultiBasis<T> mb(mp);
1073 assembler.setIntegrationElements(mb);
1074
1076 bc.setGeoMap(mp);
1077 for (auto bit = mp.bBegin(); bit != mp.bEnd(); ++bit) {
1078 bc.addCondition(*bit, condition_type::dirichlet, nullptr);
1079 }
1080
1081 space space1 = assembler.getSpace(mb, d); // 1D space!!
1082 space1.setup(bc, dirichlet::homogeneous, 0);
1083
1084 // Function for the Residual
1085 typedef std::function<gsSparseMatrix<T>(gsVector<T> const &)> Jacobian_t;
1086 typedef std::function<gsVector<T>(gsVector<T> const &)> Residual_t;
1087
1088 gsMultiPatch<T> mpSubstitute = mp;
1089 Residual_t Residual = [&assembler, &mapper, &mpSubstitute, &space1](
1090 gsVector<T> const &x) {
1091 // Convert free vector to MultiPatch
1092 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1093
1094 // Calculate geometry map
1095 geometryMap G = assembler.getMap(mpSubstitute);
1096
1097 // Compute metric matrix and Hessian of G
1098 auto metricMat = jac(G).tr() * jac(G);
1099 auto hessG = hess(G);
1100
1101 // Compute scale factor
1102 auto scale = metricMat.trace();
1103
1104 // Compute Lx
1105 auto Lx = hessG % metricMat.adj() / scale.val();
1106
1107 // Compute nonlinear system
1108 auto nonlinearSystem = frprod3(space1, Lx);
1109
1110 // Assemble the system
1111 assembler.initSystem();
1112 assembler.assemble(nonlinearSystem);
1113
1114 // Return the assembled system
1115 return assembler.rhs();
1116 };
1117
1118 Jacobian_t Jacobian;
1119 int preconditionerType = options.askInt("AAPreconditionType", 0);
1120 switch (preconditionerType) {
1121 case 1:
1122 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](gsVector<T> const &x)
1123 {
1124 // diagonal Jacobian matrix as a preconditioner
1125 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1126 geometryMap G = assembler.getMap(mpSubstitute);
1127
1128 auto jacMat = jacScaledLxDiag(space1, G);
1129
1130 assembler.initSystem();
1131 assembler.assemble(jacMat);
1132 return assembler.matrix();
1133 };
1134 break;
1135
1136 case 2:
1137 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](gsVector<T> const &x)
1138 {
1139 // diagonal-block Jacobian matrix as a preconditioner
1140 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1141 geometryMap G = assembler.getMap(mpSubstitute);
1142
1143 auto jacMat = jacScaledLxDiagBlock(space1, G);
1144
1145 assembler.initSystem();
1146 assembler.assemble(jacMat);
1147 return assembler.matrix();
1148 };
1149 break;
1150
1151 default:
1152 // analytical Jacobian matrix as a preconditioner
1153 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](gsVector<T> const &x)
1154 {
1155 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1156 geometryMap G = assembler.getMap(mpSubstitute);
1157
1158 auto jacMat = jacScaledLx(space1, G);
1159
1160 assembler.initSystem();
1161 assembler.assemble(jacMat);
1162 return assembler.matrix();
1163 };
1164 }
1165
1166 preAAParam<T> param;
1167 param.m = options.askInt("AAwindowsize", 5);
1168 param.usePreconditioning = true;
1169 param.updatePreconditionerStep = options.askInt("N_update", 5);
1170 param.epsilon = 1e-5;
1171
1172 preAApp::AndersonAcceleration<T> AASolver(param);
1173 gsVector<T> solVector = AASolver.compute(initialGuessVector,Residual, Jacobian);
1174
1175 if (options.askSwitch("needPDEH1", true))
1176 {
1177 verboseLog("\nStart parameterization improvement by H1 discrezation...",options.askInt("Verbose", 0));
1178 Residual = [&assembler, &mapper, &mpSubstitute, &space1](gsVector<T> const &x)
1179 {
1180 // H1 discretization
1181 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1182 geometryMap G = assembler.getMap(mpSubstitute);
1183
1184 auto invJacMat = jac(G).inv();
1185 auto jacU = jac(space1) * invJacMat;
1186 auto nonlinearSystem = jacU % invJacMat;
1187
1188 assembler.initSystem();
1189 assembler.assemble(nonlinearSystem);
1190 return assembler.rhs();
1191 };
1192
1193 preconditionerType = options.askInt("AAPreconditionType", 0);
1194 switch (preconditionerType) {
1195 case 1:
1196 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](gsVector<T> const &x)
1197 {
1198 // diagonal block Jacobian matrix as a preconditioner
1199 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1200 geometryMap G = assembler.getMap(mpSubstitute);
1201
1202 auto jacMat = jacScaledLxH1DiagBlock(space1, G);
1203
1204 assembler.initSystem();
1205 assembler.assemble(jacMat);
1206
1207 return assembler.matrix();
1208 };
1209 break;
1210
1211 default:
1212 Jacobian = [&assembler, &mapper, &mpSubstitute,&space1](gsVector<T> const &x)
1213 {
1214 // analytical Jacobian matrix as a preconditioner
1215 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1216 geometryMap G = assembler.getMap(mpSubstitute);
1217
1218 auto jacMat = jacScaledLxH1(space1, G);
1219
1220 assembler.initSystem();
1221 assembler.assemble(jacMat);
1222
1223 return assembler.matrix();
1224 };
1225 }
1226
1227 preAApp::AndersonAcceleration<T> AASolver2(param);
1228 solVector = AASolver2.compute(solVector, Residual, Jacobian);
1229 }
1230
1231 gsMultiPatch<T> result = mp;
1232 convertFreeVectorToMultiPatch<T>(solVector, mapper, result);
1233
1234 verboseLog("Finished!", options.askInt("Verbose", 0));
1235
1236 return result;
1237}
1238
1239namespace expr {
1240
1241template<typename E1, typename E2>
1242class frprod2_expr;
1243template<typename E1, typename E2>
1244class frprod3_expr;
1245template<class E0, class E1, class E2>
1246class ternary_expr; // ternary expression
1247template<class E>
1248class jacScaledLx_expr;
1249template<class E>
1250class jacScaledLxDiag_expr;
1251template<class E>
1252class jacScaledLxDiagBlock_expr;
1253template<class E>
1254class jacScaledLxH1_expr;
1255template<class E>
1256class jacScaledLxH1DiagBlock_expr;
1257
1258/*
1259
1260 Expression for the product of the jacob of space and matrix A. Return term
1261 is a column vector (d*n x 1). The result should be the same as the result
1262 from jac(space) % A, whereas avoids redundant multiply with zero-value components.
1263
1264 For 2D case:
1265 frprod2(space,A) = [a_{11}*dN_1/dxi_1, a_{12}*dN_1/dxi_2, a_{11}*dN_2/dxi_1,
1266 a_{12}*dN_2/dxi_2,...,a_{11}*dN_n/dxi_1, a_{12}*dN_n/dxi_2, a_{21}*dN_1/dxi_1,
1267 a_{22}*dN_1/dxi_2, a_{21}*dN_2/dxi_1, a_{22}*dN_2/dxi_2,...,a_{21}*dN_n/dxi_1,
1268 a_{22}*dN_n/dxi_2]
1269 For 3D case:
1270 frprod2(space,A) = [a_{11}*dN_1/dxi_1, a_{12}*dN_1/dxi_2, a_{13}*dN_1/dxi_3,...
1271 a_{11}*dN_2/dxi_1, a_{12}*dN_2/dxi_1, a_{13}*dN_2/dxi_3,...,
1272 a_{11}*dN_n/dxi_1, a_{12}*dN_n/dxi_1, a_{13}*dN_n/dxi_3,
1273 a_{21}*dN_1/dxi_1, a_{22}*dN_1/dxi_2, a_{23}*dN_1/dxi_3,...
1274 a_{21}*dN_2/dxi_1, a_{22}*dN_2/dxi_1, a_{23}*dN_2/dxi_3,...,
1275 a_{21}*dN_n/dxi_1, a_{22}*dN_n/dxi_1, a_{23}*dN_n/dxi_3,
1276 a_{31}*dN_1/dxi_1, a_{32}*dN_1/dxi_2, a_{33}*dN_1/dxi_3,...
1277 a_{31}*dN_2/dxi_1, a_{32}*dN_2/dxi_1, a_{33}*dN_2/dxi_3,...,
1278 a_{31}*dN_n/dxi_1, a_{32}*dN_n/dxi_1, a_{33}*dN_n/dxi_3,]
1279
1280NOTE: _u should be a space, _v should NOT be a space (fix with assert)
1281
1282 */
1283template<typename E1, typename E2>
1284class frprod3_expr : public _expr<frprod3_expr<E1, E2> > {
1285 public:
1286 typedef typename E2::Scalar Scalar;
1287 enum { ScalarValued = 0, Space = E1::Space, ColBlocks = 0 };
1288
1289 private:
1290 typename E1::Nested_t _u;
1291 typename E2::Nested_t _v;
1292
1293 mutable gsMatrix<Scalar> res, bGrads, b;
1294
1295 public:
1296
1297 frprod3_expr(_expr<E1> const &u, _expr<E2> const &v)
1298 : _u(u), _v(v) {
1299 // gsInfo << "expression is space ? "<<E1::Space <<"\n"; _u.print(gsInfo);
1300 // GISMO_ASSERT(_u.rows() == _v.rows(),
1301 // "Wrong dimensions "<<_u.rows()<<"!="<<_v.rows()<<" in % operation");
1302 // GISMO_ASSERT(_u.cols() == _v.cols(),
1303 // "Wrong dimensions "<<_u.cols()<<"!="<<_v.cols()<<" in % operation");
1304 }
1305
1306 const gsMatrix<Scalar> &eval(const index_t k) const //todo: specialize for nb==1
1307 {
1308 auto A = _v.eval(k);
1309 b = _u.eval(k);
1310// gsDebugVar(b);
1311
1312 res.noalias() = b * A;
1313 return res;
1314 }
1315
1316 index_t rows() const { return 1; }
1317 index_t cols() const { return 1; }
1318
1319 void parse(gsExprHelper<Scalar> &evList) const {
1320 _v.parse(evList);
1321// evList.add(_u);
1322// _u.data().flags |= NEED_GRAD;
1323 _u.parse(evList);
1324 }
1325
1326 const gsFeSpace<Scalar> &rowVar() const { return _u.rowVar(); }
1327 const gsFeSpace<Scalar> &colVar() const { return _v.rowVar(); }
1328
1329 void print(std::ostream &os) const {
1330 os << "(";
1331 _u.print(os);
1332 os << " % ";
1333 _v.print(os);
1334 os << ")";
1335 }
1336};
1337
1339template<typename E1, typename E2>
1340EIGEN_STRONG_INLINE
1341frprod3_expr<E1, E2> const frprod3(E1 const &u,
1342 E2 const &M) {
1343 return frprod3_expr<E1, E2>(u, M);
1344}
1345
1346template<class E0, class E1, class E2>
1347class ternary_expr : public _expr<ternary_expr<E0, E1, E2> > {
1348 typename E0::Nested_t _u;
1349 typename E1::Nested_t _v;
1350 typename E2::Nested_t _w;
1351 public:
1352 typedef typename E1::Scalar Scalar;
1353
1354 explicit ternary_expr(_expr<E0> const &u,
1355 _expr<E1> const &v,
1356 _expr<E2> const &w)
1357 :
1358 _u(u),
1359 _v(v),
1360 _w(w) {
1361 GISMO_ASSERT(E0::ScalarValued, "Condition must be scalar valued");
1362 GISMO_ASSERT((int) E1::ScalarValued == (int) E2::ScalarValued,
1363 "Both v and w must be scalar valued (or not).");
1364 GISMO_ASSERT((int) E1::ColBlocks == (int) E2::ColBlocks,
1365 "Both v and w must be colblocks (or not).");
1366 GISMO_ASSERT((int) E1::Space == (int) E2::Space,
1367 "Both v and w must be space (or not), but E1::Space = "
1368 << E1::Space << " and E2::Space = " << E2::Space);
1369 GISMO_ASSERT(_v.rows() == _w.rows(),
1370 "Rows of v and w differ. _v.rows() = " << _v.rows()
1371 << ", _w.rows() = "
1372 << _w.rows());
1373 GISMO_ASSERT(_v.cols() == _w.cols(),
1374 "Columns of v and w differ. _v.cols() = " << _v.cols()
1375 << ", _w.cols() = "
1376 << _w.cols());
1377 GISMO_ASSERT(_v.rowVar() == _w.rowVar(), "rowVar of v and w differ.");
1378 GISMO_ASSERT(_v.colVar() == _w.colVar(), "colVar of v and w differ.");
1379 }
1380 public:
1381 enum {
1382 ScalarValued = E1::ScalarValued,
1383 ColBlocks = E1::ColBlocks,
1384 Space = E1::Space
1385 }; // == E2::Space
1386
1387// const Scalar eval(const index_t k) const { return (_u.eval(k) > 0 ? _v.eval
1388// (k) : _w.eval(k)); }
1389
1390 const Temporary_t eval(const index_t k) const {
1391 return (_u.eval(k) > 0 ? _v.eval(k) : _w.eval(k));
1392 }
1393
1394 // { res = eval_impl(_u,_v,_w,k); return res;}
1395
1396 index_t rows() const { return _v.rows(); }
1397 index_t cols() const { return _v.cols(); }
1398 void parse(gsExprHelper<Scalar> &evList) const {
1399 _u.parse(evList);
1400 _v.parse(evList);
1401 _w.parse(evList);
1402 }
1403
1404 const gsFeSpace<Scalar> &rowVar() const { return _v.rowVar(); }
1405 const gsFeSpace<Scalar> &colVar() const { return _v.colVar(); }
1406
1407 void print(std::ostream &os) const { os << "ternary"; }
1408
1410// private:
1411// typename util::enable_if<E1::ScalarValued,Scalar>::type
1412// eval_impl(const E0 & u, const E1 & v, const E2 & w, const index_t k) const
1413// {
1414// return (u.eval(k) > 0 ? v.eval(k) : w.eval(k));
1415// }
1416
1417// typename util::enable_if<!ScalarValued,gsMatrix<Scalar>>::type
1418// eval_impl(const E0 & u, const E1 & v, const E2 & w, const index_t k)
1419// {
1420// return (u.eval(k).array() > 0 ? v.eval(k) : w.eval(k));
1421// }
1422
1423// private:
1424// template<class U, class V, class W> static inline
1425// typename util::enable_if<U::ScalarValued && V::ScalarValued,AutoReturn_t>::type
1426// eval_impl(const U &u, const V & v, const W & w, const index_t k)
1427// {
1428// gsMatrix<Scalar> res(1,1);
1429// // bool test = u.eval(k) > 0;
1430// res<<(u.eval(k) > 0 ? v.eval(k) : w.eval(k));
1431// return res;
1432// }
1433
1434// template<class U, class V, class W> static inline
1435// typename util::enable_if<U::ScalarValued && !V::ScalarValued,AutoReturn_t>::type
1436// eval_impl(const U &u, const V & v, const W & w, const index_t k)
1437// {
1438// return u.eval(k) > 0 ? v.eval(k) : w.eval(k);
1439// }
1440
1441// template<class U, class V, class W> static inline
1442// typename util::enable_if<!U::ScalarValued,gsMatrix<Scalar>>::type
1443// eval_impl(const U &u, const V & v, const W & w, const index_t k)
1444// {
1445// GISMO_ERROR("Something went wrong");
1446// }
1447};
1448
1449/*
1450 Expression for Jacobian matrix for PDE-based parameterization construction
1451*/
1452template<class E>
1453class jacScaledLx_expr : public _expr<jacScaledLx_expr<E> > {
1454 public:
1455 typedef typename E::Scalar Scalar;
1456
1457 private:
1458 typename E::Nested_t _u;
1459 typename gsGeometryMap<Scalar>::Nested_t _G;
1460
1461 public:
1462 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1463
1464 jacScaledLx_expr(const E &u, const gsGeometryMap<Scalar> &G) : _u(u), _G(G) {}
1465
1466 mutable gsMatrix<Scalar> res, derivGeom, deriv2Geom, derivBasis, deriv2Basis;
1467 mutable gsMatrix<Scalar> dg11dx, dg11dy, dg22dx, dg22dy, dg12dx, dg12dy;
1468 mutable gsMatrix<Scalar> commonTerm, dLxdx, dLxdy, dLydx, dLydy;
1469
1470// EIGEN_MAKE_ALIGNED_OPERATOR_NEW
1471
1472 const gsMatrix<Scalar> &eval(const index_t k) const {
1473 gsMatrix<Scalar> basis = _u.data().values[0].col(k);
1474
1475 derivBasis = _u.data().values[1].col(k).transpose();
1476 deriv2Basis = _u.data().values[2].col(k).transpose();
1477
1478 derivBasis.blockTransposeInPlace(_u.dim());
1479 deriv2Basis.blockTransposeInPlace(1 + _u.dim());
1480
1481 derivGeom = _G.data().values[1].col(k);
1482 deriv2Geom = _G.data().values[2].col(k);
1483
1484 Scalar g11 = derivGeom(0) * derivGeom(0) + derivGeom(2) * derivGeom(2);
1485 Scalar g12 = derivGeom(0) * derivGeom(1) + derivGeom(2) * derivGeom(3);
1486 Scalar g22 = derivGeom(1) * derivGeom(1) + derivGeom(3) * derivGeom(3);
1487
1488 Scalar scaleFactor = g11 + g22;
1489
1490 Scalar Lx =
1491 (g22 * deriv2Geom(0) + g11 * deriv2Geom(1) - 2.0 * g12 * deriv2Geom(2))
1492 / scaleFactor;
1493 Scalar Ly =
1494 (g22 * deriv2Geom(3) + g11 * deriv2Geom(4) - 2.0 * g12 * deriv2Geom(5))
1495 / scaleFactor;
1496
1497 dg11dx.noalias() = 2.0 * derivGeom(0) * derivBasis.row(0);
1498 dg11dy.noalias() = 2.0 * derivGeom(2) * derivBasis.row(0);
1499 dg22dx.noalias() = 2.0 * derivGeom(1) * derivBasis.row(1);
1500 dg22dy.noalias() = 2.0 * derivGeom(3) * derivBasis.row(1);
1501 dg12dx.noalias() =
1502 derivGeom(1) * derivBasis.row(0) + derivGeom(0) * derivBasis.row(1);
1503 dg12dy.noalias() =
1504 derivGeom(3) * derivBasis.row(0) + derivGeom(2) * derivBasis.row(1);
1505
1506 commonTerm.noalias() = g22 * deriv2Basis.row(0) + g11 * deriv2Basis.row(1)
1507 - 2.0 * g12 * deriv2Basis.row(2);
1508 dLxdx.noalias() = dg22dx * deriv2Geom(0) + dg11dx * deriv2Geom(1)
1509 - 2.0 * dg12dx * deriv2Geom(2) + commonTerm;
1510 dLxdy.noalias() = dg22dy * deriv2Geom(0) + dg11dy * deriv2Geom(1)
1511 - 2.0 * dg12dy * deriv2Geom(2);
1512 dLydx.noalias() = dg22dx * deriv2Geom(3) + dg11dx * deriv2Geom(4)
1513 - 2.0 * dg12dx * deriv2Geom(5);
1514 dLydy.noalias() = dg22dy * deriv2Geom(3) + dg11dy * deriv2Geom(4)
1515 - 2.0 * dg12dy * deriv2Geom(5) + commonTerm;
1516
1517 dLxdx = (dLxdx - Lx * (dg11dx + dg22dx)) / scaleFactor;
1518 dLxdy = (dLxdy - Lx * (dg11dy + dg22dy)) / scaleFactor;
1519 dLydx = (dLydx - Ly * (dg11dx + dg22dx)) / scaleFactor;
1520 dLydy = (dLydy - Ly * (dg11dy + dg22dy)) / scaleFactor;
1521
1522 const index_t A = _u.cardinality() / _u.dim(); // _u.data().actives.rows()
1523 res.resize(_u.cardinality(), _u.cardinality());
1524// res.setZero();
1525 res.topLeftCorner(A, A).noalias() = basis * dLxdx;
1526 res.topRightCorner(A, A).noalias() = basis * dLxdy;
1527 res.bottomLeftCorner(A, A).noalias() = basis * dLydx;
1528 res.bottomRightCorner(A, A).noalias() = basis * dLydy;
1529
1530 return res;
1531 }
1532
1533 index_t rows() const { return 1; }
1534
1535 index_t cols() const { return 1; }
1536
1537 void parse(gsExprHelper<Scalar> &evList) const {
1538 evList.add(_u);
1539 _u.data().flags |= NEED_VALUE | NEED_GRAD | NEED_DERIV2;
1540
1541 evList.add(_G);
1542 _G.data().flags |= NEED_DERIV | NEED_DERIV2;
1543 }
1544
1545 const gsFeSpace<Scalar> &rowVar() const { return _u.rowVar(); }
1546 const gsFeSpace<Scalar> &colVar() const { return _u.rowVar(); }
1547 // TODO: question, what do these parameters mean?
1548 index_t cardinality_impl() const { return _u.cardinality_impl(); }
1549
1550 void print(std::ostream &os) const {
1551 os << "jacScaledLx(";
1552 _u.print(os);
1553 os << ")";
1554 }
1555};
1556
1557/*
1558 Expression for Jacobian matrix (diagonal part) for PDE-based parameterization construction
1559*/
1560template<class E>
1561class jacScaledLxDiag_expr : public _expr<jacScaledLxDiag_expr<E> > {
1562 public:
1563 typedef typename E::Scalar Scalar;
1564
1565 private:
1566 typename E::Nested_t _u;
1567 typename gsGeometryMap<Scalar>::Nested_t _G;
1568
1569 public:
1570 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1571
1572 jacScaledLxDiag_expr(const E &u, const gsGeometryMap<Scalar> &G)
1573 : _u(u), _G(G) {}
1574
1575 mutable gsMatrix<Scalar> res, derivGeom, deriv2Geom, derivBasis, deriv2Basis;
1576 mutable gsMatrix<Scalar> dg11dx, dg11dy, dg22dx, dg22dy, dg12dx, dg12dy;
1577 mutable gsMatrix<Scalar> commonTerm, dLxdx, dLydy;
1578
1579// EIGEN_MAKE_ALIGNED_OPERATOR_NEW
1580
1581 const gsMatrix<Scalar> &eval(const index_t k) const {
1582 derivBasis = _u.data().values[1].col(k).transpose();
1583 deriv2Basis = _u.data().values[2].col(k).transpose();
1584
1585 derivBasis.blockTransposeInPlace(_u.dim());
1586 deriv2Basis.blockTransposeInPlace(1 + _u.dim());
1587
1588 derivGeom = _G.data().values[1].col(k);
1589 deriv2Geom = _G.data().values[2].col(k);
1590
1591 Scalar g11 = derivGeom(0) * derivGeom(0) + derivGeom(2) * derivGeom(2);
1592 Scalar g12 = derivGeom(0) * derivGeom(1) + derivGeom(2) * derivGeom(3);
1593 Scalar g22 = derivGeom(1) * derivGeom(1) + derivGeom(3) * derivGeom(3);
1594
1595 Scalar scaleFactor = g11 + g22;
1596
1597 Scalar Lx =
1598 (g22 * deriv2Geom(0) + g11 * deriv2Geom(1) - 2.0 * g12 * deriv2Geom(2))
1599 / scaleFactor;
1600 Scalar Ly =
1601 (g22 * deriv2Geom(3) + g11 * deriv2Geom(4) - 2.0 * g12 * deriv2Geom(5))
1602 / scaleFactor;
1603
1604 dg11dx.noalias() = 2.0 * derivGeom(0) * derivBasis.row(0);
1605 dg11dy.noalias() = 2.0 * derivGeom(2) * derivBasis.row(0);
1606 dg22dx.noalias() = 2.0 * derivGeom(1) * derivBasis.row(1);
1607 dg22dy.noalias() = 2.0 * derivGeom(3) * derivBasis.row(1);
1608 dg12dx.noalias() =
1609 derivGeom(1) * derivBasis.row(0) + derivGeom(0) * derivBasis.row(1);
1610 dg12dy.noalias() =
1611 derivGeom(3) * derivBasis.row(0) + derivGeom(2) * derivBasis.row(1);
1612
1613 commonTerm.noalias() = g22 * deriv2Basis.row(0) + g11 * deriv2Basis.row(1)
1614 - 2.0 * g12 * deriv2Basis.row(2);
1615 dLxdx.noalias() = dg22dx * deriv2Geom(0) + dg11dx * deriv2Geom(1)
1616 - 2.0 * dg12dx * deriv2Geom(2) + commonTerm;
1617 dLydy.noalias() = dg22dy * deriv2Geom(3) + dg11dy * deriv2Geom(4)
1618 - 2.0 * dg12dy * deriv2Geom(5) + commonTerm;
1619
1620 dLxdx = (dLxdx - Lx * (dg11dx + dg22dx)) / scaleFactor;
1621 dLydy = (dLydy - Ly * (dg11dy + dg22dy)) / scaleFactor;
1622
1623 const index_t A = _u.cardinality() / _u.dim(); // _u.data().actives.rows()
1624 res.resize(_u.cardinality(), _u.cardinality());
1625 res.setZero();
1626 res.topLeftCorner(A, A) = (_u.data().values[0].col(k).array()
1627 * dLxdx.array()).matrix().asDiagonal();
1628 res.bottomRightCorner(A, A) = (_u.data().values[0].col(k).array()
1629 * dLydy.array()).matrix().asDiagonal();
1630
1631 return res;
1632 }
1633
1634 index_t rows() const { return 1; }
1635
1636 index_t cols() const { return 1; }
1637
1638 void parse(gsExprHelper<Scalar> &evList) const {
1639 evList.add(_u);
1640 _u.data().flags |= NEED_VALUE | NEED_GRAD | NEED_DERIV2;
1641
1642 evList.add(_G);
1643 _G.data().flags |= NEED_DERIV | NEED_DERIV2;
1644 }
1645
1646 const gsFeSpace<Scalar> &rowVar() const { return _u.rowVar(); }
1647 const gsFeSpace<Scalar> &colVar() const { return _u.rowVar(); }
1648 // TODO: question, what do these parameters mean?
1649 index_t cardinality_impl() const { return _u.cardinality_impl(); }
1650
1651 void print(std::ostream &os) const {
1652 os << "jacScaledLxDiag(";
1653 _u.print(os);
1654 os << ")";
1655 }
1656};
1657
1658/*
1659 Expression for Jacobian matrix (diagonal block) for PDE-based parameterization construction
1660*/
1661template<class E>
1662class jacScaledLxDiagBlock_expr : public _expr<jacScaledLxDiagBlock_expr<E> > {
1663 public:
1664 typedef typename E::Scalar Scalar;
1665
1666 private:
1667 typename E::Nested_t _u;
1668 typename gsGeometryMap<Scalar>::Nested_t _G;
1669
1670 public:
1671 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1672
1673 jacScaledLxDiagBlock_expr(const E &u, const gsGeometryMap<Scalar> &G)
1674 : _u(u), _G(G) {}
1675
1676 mutable gsMatrix<Scalar> res, basis, derivGeom, deriv2Geom, derivBasis,
1677 deriv2Basis;
1678 mutable gsMatrix<Scalar> dg11dx, dg11dy, dg22dx, dg22dy, dg12dx, dg12dy;
1679 mutable gsMatrix<Scalar> commonTerm, dLxdx, dLydy;
1680
1681// EIGEN_MAKE_ALIGNED_OPERATOR_NEW
1682
1683 const gsMatrix<Scalar> &eval(const index_t k) const {
1684// basis = _u.data().values[0].col(k);
1685
1686 derivBasis = _u.data().values[1].col(k).transpose();
1687 deriv2Basis = _u.data().values[2].col(k).transpose();
1688
1689 derivBasis.blockTransposeInPlace(_u.dim());
1690 deriv2Basis.blockTransposeInPlace(1 + _u.dim());
1691
1692 derivGeom = _G.data().values[1].col(k);
1693 deriv2Geom = _G.data().values[2].col(k);
1694
1695 Scalar g11 = derivGeom(0) * derivGeom(0) + derivGeom(2) * derivGeom(2);
1696 Scalar g12 = derivGeom(0) * derivGeom(1) + derivGeom(2) * derivGeom(3);
1697 Scalar g22 = derivGeom(1) * derivGeom(1) + derivGeom(3) * derivGeom(3);
1698 Scalar scaleFactor = g11 + g22;
1699
1700 Scalar Lx =
1701 (g22 * deriv2Geom(0) + g11 * deriv2Geom(1) - 2.0 * g12 * deriv2Geom(2))
1702 / scaleFactor;
1703 Scalar Ly =
1704 (g22 * deriv2Geom(3) + g11 * deriv2Geom(4) - 2.0 * g12 * deriv2Geom(5))
1705 / scaleFactor;
1706
1707 dg11dx.noalias() = 2.0 * derivGeom(0) * derivBasis.row(0);
1708 dg11dy.noalias() = 2.0 * derivGeom(2) * derivBasis.row(0);
1709 dg22dx.noalias() = 2.0 * derivGeom(1) * derivBasis.row(1);
1710 dg22dy.noalias() = 2.0 * derivGeom(3) * derivBasis.row(1);
1711 dg12dx.noalias() =
1712 derivGeom(1) * derivBasis.row(0) + derivGeom(0) * derivBasis.row(1);
1713 dg12dy.noalias() =
1714 derivGeom(3) * derivBasis.row(0) + derivGeom(2) * derivBasis.row(1);
1715
1716 commonTerm.noalias() = g22 * deriv2Basis.row(0) + g11 * deriv2Basis.row(1)
1717 - 2.0 * g12 * deriv2Basis.row(2);
1718 dLxdx.noalias() = dg22dx * deriv2Geom(0) + dg11dx * deriv2Geom(1)
1719 - 2.0 * dg12dx * deriv2Geom(2) + commonTerm;
1720 dLydy.noalias() = dg22dy * deriv2Geom(3) + dg11dy * deriv2Geom(4)
1721 - 2.0 * dg12dy * deriv2Geom(5) + commonTerm;
1722
1723 dLxdx = (dLxdx - Lx * (dg11dx + dg22dx)) / scaleFactor;
1724 dLydy = (dLydy - Ly * (dg11dy + dg22dy)) / scaleFactor;
1725
1726 const index_t A = _u.cardinality() / _u.dim(); // _u.data().actives.rows()
1727 res.resize(_u.cardinality(), _u.cardinality());
1728 res.setZero();
1729 res.topLeftCorner(A, A).noalias() =
1730 _u.data().values[0].col(k) * dLxdx;
1731 res.bottomRightCorner(A, A).noalias() =
1732 _u.data().values[0].col(k) * dLydy;
1733
1734 return res;
1735 }
1736
1737 index_t rows() const { return 1; }
1738
1739 index_t cols() const { return 1; }
1740
1741 void parse(gsExprHelper<Scalar> &evList) const {
1742 evList.add(_u);
1743 _u.data().flags |= NEED_VALUE | NEED_GRAD | NEED_DERIV2;
1744
1745 evList.add(_G);
1746 _G.data().flags |= NEED_DERIV | NEED_DERIV2;
1747 }
1748
1749 const gsFeSpace<Scalar> &rowVar() const { return _u.rowVar(); }
1750 const gsFeSpace<Scalar> &colVar() const { return _u.rowVar(); }
1751 // TODO: question, what do these parameters mean?
1752 index_t cardinality_impl() const { return _u.cardinality_impl(); }
1753
1754 void print(std::ostream &os) const {
1755 os << "jacScaledLxDiagBlock(";
1756 _u.print(os);
1757 os << ")";
1758 }
1759};
1760
1761/*
1762 Expression for Jacobian matrix (in H1 space) for PDE-based parameterization construction
1763*/
1764template<class E>
1765class jacScaledLxH1_expr : public _expr<jacScaledLxH1_expr<E> > {
1766 public:
1767 typedef typename E::Scalar Scalar;
1768
1769 private:
1770 typename E::Nested_t _u;
1771 typename gsGeometryMap<Scalar>::Nested_t _G;
1772
1773 public:
1774 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1775
1776 jacScaledLxH1_expr(const E &u, const gsGeometryMap<Scalar> &G)
1777 : _u(u), _G(G) {}
1778
1779 mutable gsMatrix<Scalar> res, derivBasis;
1780
1781// EIGEN_MAKE_ALIGNED_OPERATOR_NEW
1782
1783 const gsMatrix<Scalar> &eval(const index_t k) const {
1784 gsMatrix<Scalar> jacMat = _G.data().values[1].reshapeCol(k,
1785 _G.data().dim.first,
1786 _G.data().dim.second);
1787 gsMatrix<Scalar> invJacMat = jacMat.inverse();
1788
1789 derivBasis = _u.data().values[1].col(k).transpose();
1790 derivBasis.blockTransposeInPlace(_u.dim());
1791 gsMatrix<Scalar> invJacHatG;
1792 invJacHatG.noalias() = invJacMat * derivBasis;
1793
1794 const index_t N = _u.cardinality() / _u.dim(); // _u.data().actives.rows()
1795 gsMatrix<Scalar> jacdLxdx(N, N);
1796 gsMatrix<Scalar> jacdLxdy(N, N);
1797 gsMatrix<Scalar> jacdLydx(N, N);
1798 gsMatrix<Scalar> jacdLydy(N, N);
1799
1800 gsMatrix<> temp(2, N);
1801 for (auto i = 0; i < N; ++i) {
1802 // for x-direction
1803 temp.row(0).noalias() = invJacHatG(0, i) * invJacHatG.row(0);
1804 temp.row(1).noalias() = invJacHatG(1, i) * invJacHatG.row(0);
1805 gsMatrix<Scalar> dinvJacdx(2, 2);
1806 dinvJacdx.row(0).noalias() = invJacHatG(0, i) * invJacMat.row(0);
1807 dinvJacdx.row(1).noalias() = invJacHatG(1, i) * invJacMat.row(0);
1808 jacdLxdx.col(i).noalias() = -temp.transpose() * invJacMat.col(0)
1809 - invJacHatG.transpose() * dinvJacdx.col(0);
1810 jacdLydx.col(i).noalias() = -temp.transpose() * invJacMat.col(1)
1811 - invJacHatG.transpose() * dinvJacdx.col(1);
1812
1813 // for y-direction
1814 temp.row(0).noalias() = invJacHatG(0, i) * invJacHatG.row(1);
1815 temp.row(1).noalias() = invJacHatG(1, i) * invJacHatG.row(1);
1816 gsMatrix<Scalar> dinvJacdy(2, 2);
1817 dinvJacdy.row(0).noalias() = invJacHatG(0, i) * invJacMat.row(1);
1818 dinvJacdy.row(1).noalias() = invJacHatG(1, i) * invJacMat.row(1);
1819 jacdLxdy.col(i).noalias() = -temp.transpose() * invJacMat.col(0)
1820 - invJacHatG.transpose() * dinvJacdy.col(0);
1821 jacdLydy.col(i).noalias() = -temp.transpose() * invJacMat.col(1)
1822 - invJacHatG.transpose() * dinvJacdy.col(1);
1823 }
1824
1825 res.resize(_u.cardinality(), _u.cardinality());
1826 res.setZero();
1827 res.topLeftCorner(N, N).noalias() = jacdLxdx;
1828 res.topRightCorner(N, N).noalias() = jacdLxdy;
1829 res.bottomLeftCorner(N, N).noalias() = jacdLydx;
1830 res.bottomRightCorner(N, N).noalias() = jacdLydy;
1831
1832 return res;
1833 }
1834
1835 index_t rows() const { return 1; }
1836
1837 index_t cols() const { return 1; }
1838
1839 void parse(gsExprHelper<Scalar> &evList) const {
1840 evList.add(_u);
1841 _u.data().flags |= NEED_GRAD;
1842
1843 evList.add(_G);
1844 _G.data().flags |= NEED_DERIV;
1845 }
1846
1847 const gsFeSpace<Scalar> &rowVar() const { return _u.rowVar(); }
1848 const gsFeSpace<Scalar> &colVar() const { return _u.rowVar(); }
1849 // TODO: question, what do these parameters mean?
1850 index_t cardinality_impl() const { return _u.cardinality_impl(); }
1851
1852 void print(std::ostream &os) const {
1853 os << "jacScaledLx(";
1854 _u.print(os);
1855 os << ")";
1856 }
1857};
1858
1859/*
1860 Expression for Jacobian matrix (in H1 space) for PDE-based parameterization construction
1861*/
1862template<class E>
1863class jacScaledLxH1DiagBlock_expr
1864 : public _expr<jacScaledLxH1DiagBlock_expr<E> > {
1865 public:
1866 typedef typename E::Scalar Scalar;
1867
1868 private:
1869 typename E::Nested_t _u;
1870 typename gsGeometryMap<Scalar>::Nested_t _G;
1871
1872 public:
1873 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1874
1875 jacScaledLxH1DiagBlock_expr(const E &u, const gsGeometryMap<Scalar> &G)
1876 : _u(u), _G(G) {}
1877
1878 mutable gsMatrix<Scalar> res, derivBasis;
1879
1880// EIGEN_MAKE_ALIGNED_OPERATOR_NEW
1881
1882 const gsMatrix<Scalar> &eval(const index_t k) const {
1883 gsMatrix<Scalar> jacMat = _G.data().values[1].reshapeCol(k,
1884 _G.data().dim.first,
1885 _G.data().dim.second);
1886 gsMatrix<Scalar> invJacMat = jacMat.inverse();
1887
1888 derivBasis = _u.data().values[1].col(k).transpose();
1889 derivBasis.blockTransposeInPlace(_u.dim());
1890 gsMatrix<Scalar> invJacHatG = invJacMat * derivBasis;
1891
1892 const index_t N = _u.cardinality() / _u.dim(); // _u.data().actives.rows()
1893 gsMatrix<Scalar> jacdLxdx(N, N);
1894// gsMatrix<Scalar> jacdLxdy(N,N);
1895// gsMatrix<Scalar> jacdLydx(N,N);
1896 gsMatrix<Scalar> jacdLydy(N, N);
1897
1898 gsMatrix<> temp(2, N);
1899 for (auto i = 0; i < N; ++i) {
1900 // for x-direction
1901 temp.row(0) = invJacHatG(0, i) * invJacHatG.row(0);
1902 temp.row(1) = invJacHatG(1, i) * invJacHatG.row(0);
1903 gsMatrix<Scalar> dinvJacdx(2, 2);
1904 dinvJacdx.row(0) = invJacHatG(0, i) * invJacMat.row(0);
1905 dinvJacdx.row(1) = invJacHatG(1, i) * invJacMat.row(0);
1906 jacdLxdx.col(i) = -temp.transpose() * invJacMat.col(0)
1907 - invJacHatG.transpose() * dinvJacdx.col(0);
1908// jacdLydx.col(i) = -temp.transpose()*invJacMat.col(1) - invJacHatG.transpose()*dinvJacdx.col(1);
1909
1910 // for y-direction
1911 temp.row(0) = invJacHatG(0, i) * invJacHatG.row(1);
1912 temp.row(1) = invJacHatG(1, i) * invJacHatG.row(1);
1913 gsMatrix<Scalar> dinvJacdy(2, 2);
1914 dinvJacdy.row(0) = invJacHatG(0, i) * invJacMat.row(1);
1915 dinvJacdy.row(1) = invJacHatG(1, i) * invJacMat.row(1);
1916// jacdLxdy.col(i) = -temp.transpose()*invJacMat.col(0) - invJacHatG.transpose()*dinvJacdy.col(0);
1917 jacdLydy.col(i) = -temp.transpose() * invJacMat.col(1)
1918 - invJacHatG.transpose() * dinvJacdy.col(1);
1919 }
1920
1921 res.resize(_u.cardinality(), _u.cardinality());
1922 res.setZero();
1923 res.topLeftCorner(N, N) = jacdLxdx;
1924// res.topRightCorner(N,N) = jacdLxdy;
1925// res.bottomLeftCorner(N,N) = jacdLydx;
1926 res.bottomRightCorner(N, N) = jacdLydy;
1927
1928 return res;
1929 }
1930
1931 index_t rows() const { return 1; }
1932
1933 index_t cols() const { return 1; }
1934
1935 void parse(gsExprHelper<Scalar> &evList) const {
1936 evList.add(_u);
1937 _u.data().flags |= NEED_GRAD;
1938
1939 evList.add(_G);
1940 _G.data().flags |= NEED_DERIV;
1941 }
1942
1943 const gsFeSpace<Scalar> &rowVar() const { return _u.rowVar(); }
1944 const gsFeSpace<Scalar> &colVar() const { return _u.rowVar(); }
1945 // TODO: question, what do these parameters mean?
1946 index_t cardinality_impl() const { return _u.cardinality_impl(); }
1947
1948 void print(std::ostream &os) const {
1949 os << "jacScaledLx(";
1950 _u.print(os);
1951 os << ")";
1952 }
1953};
1954
1955/*
1956
1957 Expression for the product of the jacob of space and matrix A. Return term
1958 is a column vector (d*n x 1). The result should be the same as the result
1959 from jac(space) % A, whereas avoids redundant multiply with zero-value components.
1960
1961 For 2D case:
1962 frprod2(space,A) = [a_{11}*dN_1/dxi_1, a_{12}*dN_1/dxi_2, a_{11}*dN_2/dxi_1,
1963 a_{12}*dN_2/dxi_2,...,a_{11}*dN_n/dxi_1, a_{12}*dN_n/dxi_2, a_{21}*dN_1/dxi_1,
1964 a_{22}*dN_1/dxi_2, a_{21}*dN_2/dxi_1, a_{22}*dN_2/dxi_2,...,a_{21}*dN_n/dxi_1,
1965 a_{22}*dN_n/dxi_2]
1966 For 3D case:
1967 frprod2(space,A) = [a_{11}*dN_1/dxi_1, a_{12}*dN_1/dxi_2, a_{13}*dN_1/dxi_3,...
1968 a_{11}*dN_2/dxi_1, a_{12}*dN_2/dxi_1, a_{13}*dN_2/dxi_3,...,
1969 a_{11}*dN_n/dxi_1, a_{12}*dN_n/dxi_1, a_{13}*dN_n/dxi_3,
1970 a_{21}*dN_1/dxi_1, a_{22}*dN_1/dxi_2, a_{23}*dN_1/dxi_3,...
1971 a_{21}*dN_2/dxi_1, a_{22}*dN_2/dxi_1, a_{23}*dN_2/dxi_3,...,
1972 a_{21}*dN_n/dxi_1, a_{22}*dN_n/dxi_1, a_{23}*dN_n/dxi_3,
1973 a_{31}*dN_1/dxi_1, a_{32}*dN_1/dxi_2, a_{33}*dN_1/dxi_3,...
1974 a_{31}*dN_2/dxi_1, a_{32}*dN_2/dxi_1, a_{33}*dN_2/dxi_3,...,
1975 a_{31}*dN_n/dxi_1, a_{32}*dN_n/dxi_1, a_{33}*dN_n/dxi_3,]
1976
1977NOTE: _u should be a space, _v should NOT be a space (fix with assert)
1978
1979 */
1980template<typename E1, typename E2>
1981class frprod2_expr : public _expr<frprod2_expr<E1, E2> > {
1982 public:
1983 typedef typename E2::Scalar Scalar;
1984 enum { ScalarValued = 0, Space = E1::Space, ColBlocks = 0 };
1985
1986 private:
1987 typename E1::Nested_t _u;
1988 typename E2::Nested_t _v;
1989
1990 mutable gsMatrix<Scalar> res, bGrads;
1991
1992 public:
1993
1994 frprod2_expr(_expr<E1> const &u, _expr<E2> const &v)
1995 : _u(u), _v(v) {
1996 // gsInfo << "expression is space ? "<<E1::Space <<"\n"; _u.print(gsInfo);
1997 // GISMO_ASSERT(_u.rows() == _v.rows(),
1998 // "Wrong dimensions "<<_u.rows()<<"!="<<_v.rows()<<" in % operation");
1999 // GISMO_ASSERT(_u.cols() == _v.cols(),
2000 // "Wrong dimensions "<<_u.cols()<<"!="<<_v.cols()<<" in % operation");
2001 }
2002
2003 const gsMatrix<Scalar> &eval(const index_t k) const //todo: specialize for nb==1
2004 {
2005 auto A = _v.eval(k);
2006 bGrads = _u.data().values[1].col(k).transpose();
2007 bGrads.blockTransposeInPlace(_u.dim());
2008 res.noalias() = A * bGrads;
2009 res.transposeInPlace();
2010 res.resize(1, _u.cardinality());
2011// res.reshaped(1,_u.cardinality());
2012 res.transposeInPlace();
2013 return res;
2014 }
2015
2016 index_t rows() const { return 1; }
2017 index_t cols() const { return 1; }
2018
2019 void parse(gsExprHelper<Scalar> &evList) const {
2020 _v.parse(evList);
2021 evList.add(_u);
2022 _u.data().flags |= NEED_GRAD;
2023 }
2024
2025 const gsFeSpace<Scalar> &rowVar() const { return _u.rowVar(); }
2026 const gsFeSpace<Scalar> &colVar() const { return _v.rowVar(); }
2027
2028 void print(std::ostream &os) const {
2029 os << "(";
2030 _u.print(os);
2031 os << " % ";
2032 _v.print(os);
2033 os << ")";
2034 }
2035};
2036
2038template<class E>
2039EIGEN_STRONG_INLINE
2040jacScaledLx_expr<E> jacScaledLx(const E &u,
2041 const gsGeometryMap<typename E::Scalar> &G) {
2042 return jacScaledLx_expr<E>(u, G);
2043}
2044
2046template<class E>
2047EIGEN_STRONG_INLINE
2048jacScaledLxDiag_expr<E> jacScaledLxDiag(const E &u,
2049 const gsGeometryMap<typename E::Scalar> &G) {
2050 return jacScaledLxDiag_expr<E>(u, G);
2051}
2052
2054template<class E>
2055EIGEN_STRONG_INLINE
2056jacScaledLxDiagBlock_expr<E> jacScaledLxDiagBlock(const E &u,
2057 const gsGeometryMap<typename E::Scalar> &G) {
2058 return jacScaledLxDiagBlock_expr<E>(u, G);
2059}
2060
2062template<class E>
2063EIGEN_STRONG_INLINE
2064jacScaledLxH1_expr<E> jacScaledLxH1(const E &u,
2065 const gsGeometryMap<typename E::Scalar> &G) {
2066 return jacScaledLxH1_expr<E>(u, G);
2067}
2068
2070template<class E>
2071EIGEN_STRONG_INLINE
2072jacScaledLxH1DiagBlock_expr<E> jacScaledLxH1DiagBlock(const E &u,
2073 const gsGeometryMap<
2074 typename E::Scalar> &G) {
2075 return jacScaledLxH1DiagBlock_expr<E>(u, G);
2076}
2077
2078// Frobenious product (also known as double dot product) operator for expressions
2079template<typename E1, typename E2>
2080EIGEN_STRONG_INLINE
2081frprod2_expr<E1, E2> const frprod2(E1 const &u,
2082 E2 const &M) {
2083 return frprod2_expr<E1, E2>(u, M);
2084}
2085
2087template<class E0, class E1, class E2>
2088EIGEN_STRONG_INLINE
2089ternary_expr<E0, E1, E2> ternary(const E0 &u,
2090 const E1 &v,
2091 const E2 &w) {
2092 return ternary_expr<E0, E1, E2>(u, v, w);
2093}
2094} // namespace expr
2095
2096}// namespace gismo
Definition gsExpressions.h:973
static gsMultiPatch< T > computePenaltyPatch(const gsMultiPatch< T > &mp, const gsDofMapper &mapper, const gsOptionList &options)
Penalty function-based method (1)
static gsMultiPatch< T > computePenaltyPatch2(const gsMultiPatch< T > &mp, const gsDofMapper &mapper, const gsOptionList &options)
Penalty function-based method (2)
static gsMultiPatch< T > compute(const gsMultiPatch< T > &mp, const gsDofMapper &mapper, const gsOptionList &options)
construct analysis-suitable parameterization
Definition gsBarrierCore.hpp:131
static T computeAreaInterior(const gsMultiPatch< T > &mp)
Computes the area of a multipatch representing a full domain.
Definition gsBarrierCore.hpp:96
static T computeArea(const gsMultiPatch< T > &mp)
Compute the area of the computational domain.
Definition gsBarrierCore.hpp:90
static gsMultiPatch< T > computeVHPatch(const gsMultiPatch< T > &mp, const gsDofMapper &mapper, const gsOptionList &options)
variational harmonic method
static T computeAreaBoundary(const gsMultiPatch< T > &mp)
Computes the area of a multipatch representing boundary curves.
Definition gsBarrierCore.hpp:116
static gsMultiPatch< T > computePDEPatch(const gsMultiPatch< T > &mp, const gsDofMapper &mapper, const gsOptionList &options)
PDE-based methods, including H2 and H1.
Definition gsBarrierCore.hpp:1048
static gsMultiPatch< T > computeBarrierPatch(const gsMultiPatch< T > &mp, const gsDofMapper &mapper, const gsOptionList &options)
Barrier function-based method.
static gsOptionList defaultOptions()
Default options.
Definition gsBarrierCore.hpp:26
Class containing a set of boundary conditions.
Definition gsBoundaryConditions.h:342
void addCondition(int p, boxSide s, condition_type::type t, gsFunctionSet< T > *f, short_t unknown=0, bool parametric=false, int comp=-1)
Adds another boundary condition.
Definition gsBoundaryConditions.h:650
void setGeoMap(const gsFunctionSet< T > &gm)
Set the geometry map to evaluate boundary conditions.
Definition gsBoundaryConditions.h:916
Maintains a mapping from patch-local dofs to global dof indices and allows the elimination of individ...
Definition gsDofMapper.h:69
Definition gsExprAssembler.h:33
Generic evaluator of isogeometric expressions.
Definition gsExprEvaluator.h:39
void setIntegrationElements(const gsMultiBasis< T > &mesh)
Sets the domain of integration.
Definition gsExprEvaluator.h:110
T integral(const expr::_expr< E > &expr)
Calculates the integral of the expression expr on the whole integration domain.
Definition gsExprEvaluator.h:154
geometryMap getMap(const gsMultiPatch< T > &mp)
Registers mp as an isogeometric geometry map and return a handle to it.
Definition gsExprEvaluator.h:116
Holds a set of patch-wise bases and their topology information.
Definition gsMultiBasis.h:37
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition gsMultiPatch.h:100
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:201
bool askSwitch(const std::string &label, const bool &value=false) const
Reads value for option label from options.
Definition gsOptionList.cpp:128
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:211
index_t askInt(const std::string &label, const index_t &value=0) const
Reads value for option label from options.
Definition gsOptionList.cpp:117
void addSwitch(const std::string &label, const std::string &desc, const bool &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:235
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Definition gsBarrierCore.h:504
int m
Definition gsBarrierCore.h:513
int updatePreconditionerStep
Definition gsBarrierCore.h:542
Scalar epsilon
Definition gsBarrierCore.h:521
bool usePreconditioning
Definition gsBarrierCore.h:552
Anderson acceleration solver and its (preconditioned) variants.
Definition gsBarrierCore.h:607
const gsVector< Scalar > & compute(const gsVector< Scalar > &u0, const Residual_t &F, const Jacobian_t &Jacobian)
perform Anderson acceleration iteration
Definition gsBarrierCore.h:618
gsVector< T > convertMultiPatchToFreeVector(const gsMultiPatch< T > &mp, const gsDofMapper &mapper)
Computes a patch parametrization given a set of boundary geometries. Parametrization is not guarantee...
Definition gsBarrierCore.h:50
#define index_t
Definition gsConfig.h:32
Provides declaration of ConstantFunction class.
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
#define GISMO_ERROR(message)
Definition gsDebug.h:118
#define gsWarn
Definition gsDebug.h:50
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
void derivBasis()
Definition gsBSplineAlgorithms.h:199
EIGEN_STRONG_INLINE ternary_expr< E0, E1, E2 > ternary(const E0 &u, const E1 &v, const E2 &w)
Ternary ternary_expr.
Definition gsBarrierCore.hpp:2089
EIGEN_STRONG_INLINE jacScaledLxH1_expr< E > jacScaledLxH1(const E &u, const gsGeometryMap< typename E::Scalar > &G)
jacobian matrix of scaled Lx (in H1 space) for PDE-based parameterization construction
Definition gsBarrierCore.hpp:2064
EIGEN_STRONG_INLINE jacScaledLxH1DiagBlock_expr< E > jacScaledLxH1DiagBlock(const E &u, const gsGeometryMap< typename E::Scalar > &G)
jacobian matrix of scaled Lx (in H1 space) for PDE-based parameterization construction
Definition gsBarrierCore.hpp:2072
EIGEN_STRONG_INLINE jacScaledLxDiag_expr< E > jacScaledLxDiag(const E &u, const gsGeometryMap< typename E::Scalar > &G)
diagonal part of jacobian matrix of scaled Lx for PDE-based parameterization construction
Definition gsBarrierCore.hpp:2048
EIGEN_STRONG_INLINE frprod3_expr< E1, E2 > const frprod3(E1 const &u, E2 const &M)
Frobenious product (also known as double dot product) operator for expressions.
Definition gsBarrierCore.hpp:1341
EIGEN_STRONG_INLINE jacScaledLx_expr< E > jacScaledLx(const E &u, const gsGeometryMap< typename E::Scalar > &G)
jacobian matrix of scaled Lx for PDE-based parameterization construction
Definition gsBarrierCore.hpp:2040
EIGEN_STRONG_INLINE jac_expr< E > jac(const symbol_expr< E > &u)
The Jacobian matrix of a FE variable.
Definition gsExpressions.h:4528
EIGEN_STRONG_INLINE jacScaledLxDiagBlock_expr< E > jacScaledLxDiagBlock(const E &u, const gsGeometryMap< typename E::Scalar > &G)
diagonal block part of jacobian matrix of scaled Lx for PDE-based parameterization construction
Definition gsBarrierCore.hpp:2056
The G+Smo namespace, containing all definitions for the library.
void verboseLog(const std::string &message, const index_t &verbose)
helper function to set optimizer options
Definition gsBarrierCore.h:101
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ NEED_DERIV2
Second derivatives.
Definition gsForwardDeclarations.h:80
@ NEED_DERIV
Gradient of the object.
Definition gsForwardDeclarations.h:73
void convertFreeVectorToMultiPatch(const gsVector< T > &gsFreeVec, const gsDofMapper &mapper, gsMultiPatch< T > &mp)
Convert free control points from a vector into a multi-patch.
Definition gsBarrierCore.h:70
@ dirichlet
Dirichlet type.
Definition gsBoundaryConditions.h:31