25template<
short_t d,
typename T>
30 options.
addInt(
"Verbose",
"Print output level", 0);
33 options.
addInt(
"InitialMethod",
"Initialization method", 0);
36 options.
addInt(
"ParamMethod",
"Parameterization method", 0);
39 options.
addReal(
"ff_Delta",
"Delta for foldover-free optimization", 0.05);
40 options.
addInt(
"ff_MaxIterations",
41 "Max iterations for quality improvement",
43 options.
addReal(
"ff_MinGradientLength",
44 "Min gradient length for foldover-free optimization",
46 options.
addReal(
"ff_MinStepLength",
47 "Min step length for foldover-free optimization",
51 options.
addInt(
"qi_MaxIterations",
52 "Max iterations for quality improvement",
54 options.
addReal(
"qi_MinGradientLength",
55 "Min gradient length for quality improvement",
57 options.
addReal(
"qi_MinStepLength",
58 "Min step length for quality improvement",
63 "Quadrature points: quA*deg + quB; For patchRule: Order of target space",
66 "Quadrature points: quA*deg + quB; For patchRule: Regularity of target space",
69 "Quadrature rule [1:GaussLegendre,2:GaussLobatto,3:PatchRule]",
71 options.
addInt(
"overInt",
"Apply over-integration?", 0);
74 options.
addInt(
"AAPreconditionType",
"Preconditioner type for AA", 0);
77 options.
addInt(
"N_update",
"update preconditioner every N_update steps", 10);
80 options.
addInt(
"AAwindowsize",
"window size for preconditioned AA solver", 5);
83 options.
addSwitch(
"needPDEH1",
"improve quality by H1 discrezation?",
true);
89template<
short_t d,
typename T>
91 return computeAreaInterior(mp);
95template<
short_t d,
typename T>
107 geometryMap geometry = evaluator.
getMap(multiPatch);
110 T area = evaluator.
integral(jac(geometry).det());
115template<
short_t d,
typename T>
121enum class ParamMethod {
126 VariationalHarmonicPatch,
130template<
short_t d,
typename T>
137 method =
static_cast<ParamMethod
>(options.
askInt(
"ParamMethod", 0));
140 case ParamMethod::PenaltyPatch: {
141#ifdef gsHLBFGS_ENABLED
142 result = computePenaltyPatch(mp, mapper, options);
144 GISMO_ERROR(
"PenaltyPatch not available without gsHLBFGS. Please "
145 "compile with gsHLBGFS enabled");
149 case ParamMethod::PenaltyPatch2: {
150#ifdef gsHLBFGS_ENABLED
151 result = computePenaltyPatch2(mp, mapper, options);
153 GISMO_ERROR(
"PenaltyPatch2 not available without gsHLBFGS. Please "
154 "compile with gsHLBGFS enabled");
158 case ParamMethod::BarrierPatch: {
159#ifdef gsHLBFGS_ENABLED
160 result = computeBarrierPatch(mp, mapper, options);
162 GISMO_ERROR(
"BarrierPatch not available without gsHLBFGS. Please "
163 "compile with gsHLBGFS enabled");
167 case ParamMethod::VariationalHarmonicPatch: {
168#ifdef gsHLBFGS_ENABLED
169 result = computeVHPatch(mp, mapper, options);
171 GISMO_ERROR(
"VHPatch not available without gsHLBFGS. Please "
172 "compile with gsHLBGFS enabled");
176 case ParamMethod::PDEPatch:
178 if (method != ParamMethod::PDEPatch) {
179 gsWarn <<
"Invalid ParamMethod value. Defaulting to PDEPatch.\n";
181 result = computePDEPatch(mp, mapper, options);
188#ifdef gsHLBFGS_ENABLED
190template<
short_t d,
typename T>
196 verboseLog(
"Start variational harmonic parameterization construction...",
197 options.
askInt(
"Verbose", 0));
202 gsObjVHPt<d, T> objVHPt(mp, mapper);
203 objVHPt.applyOptions(options);
205 gsHLBFGS<T> optimizer(&objVHPt);
206 setOptimizerOptions<T>(optimizer, options);
208 optimizer.solve(initialGuessVector);
210 convertFreeVectorToMultiPatch<T>(optimizer.currentDesign(), mapper, result);
217template<
short_t d,
typename T>
220 const gsDofMapper &mapper,
221 const gsOptionList &options) {
223 verboseLog(
"Start penalty function-based parameterization construction...",
224 options.askInt(
"Verbose", 0));
227 T scaledArea = computeAreaInterior(mp);
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);
240 gsHLBFGS<T> optimizer(&objPenaltyPt);
241 setOptimizerOptions<T>(optimizer, options);
243 optimizer.solve(initialGuessVector);
244 gsMultiPatch<T> result = mp;
245 convertFreeVectorToMultiPatch<T>(optimizer.currentDesign(), mapper, result);
247 verboseLog(
"Finished!", options.askInt(
"Verbose", 0));
252template<
short_t d,
typename T>
255 const gsDofMapper &mapper,
256 const gsOptionList &options) {
258 verboseLog(
"Penalty function-based (2) parameterization construction ...",
259 options.askInt(
"Verbose", 0));
262 T scaledArea = computeAreaInterior(mp);
267 gsObjPenaltyPt2<d, T> objPenaltyPt(mp, mapper);
268 objPenaltyPt.setEps(options.askReal(
"penaltyFactor", 0.01) * scaledArea);
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);
277 gsHLBFGS<T> optimizer(&objPenaltyPt);
278 setOptimizerOptions<T>(optimizer, options);
280 optimizer.solve(initialGuessVector);
281 gsMultiPatch<T> result = mp;
282 convertFreeVectorToMultiPatch<T>(optimizer.currentDesign(), mapper, result);
284 verboseLog(
"Finished!", options.askInt(
"Verbose", 0));
289template<
short_t d,
typename T>
292 const gsDofMapper &mapper,
293 const gsOptionList &options) {
296 T scaledArea = computeAreaInterior(mp);
302 foldoverElimination(mp, mapper, initialGuessVector, scaledArea, options);
305 qualityImprovement(mp, mapper, initialGuessVector, scaledArea, options);
308 gsMultiPatch<T> result = mp;
314template<
short_t d,
typename T>
315void gsBarrierCore<d, T>::foldoverElimination(
const gsMultiPatch<T> &mp,
316 const gsDofMapper &mapper,
317 gsVector<T> &initialGuessVector,
319 const gsOptionList &options) {
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);
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",
334 optFoldoverFree.options().setReal(
"MinStepLength",
335 options.askReal(
"ff_MinStepLength", 1e-12));
336 optFoldoverFree.options().setInt(
"Verbose", options.askInt(
"Verbose", 0));
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;
341 objFoldoverFree.setDelta(delta);
343 optFoldoverFree.solve(initialGuessVector);
345 Efoldover = optFoldoverFree.objective();
346 initialGuessVector = optFoldoverFree.currentDesign();
348 if (Efoldover <= EPSILON) {
break; }
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.");
359template<
short_t d,
typename T>
360void gsBarrierCore<d, T>::qualityImprovement(
const gsMultiPatch<T> &mp,
361 const gsDofMapper &mapper,
362 gsVector<T> &initialGuessVector,
364 const gsOptionList &options) {
365 verboseLog(
"Start parameterization quality improvement step...",
366 options.askInt(
"Verbose", 0));
367 gsObjQualityImprovePt<d, T> objQualityImprovePt(mp, mapper);
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);
376 gsHLBFGS<T> optQualityImprovePt(&objQualityImprovePt);
377 setOptimizerOptions<T>(optQualityImprovePt, options);
379 optQualityImprovePt.solve(initialGuessVector);
380 initialGuessVector = optQualityImprovePt.currentDesign();
383template<
short_t d,
typename T>
384gsObjFoldoverFree<d, T>::gsObjFoldoverFree(
const gsMultiPatch<T> &patches,
388 m_mapper(std::move(mapper)),
391 m_assembler.setIntegrationElements(m_mb);
392 m_evaluator = gsExprEvaluator<T>(m_assembler);
395template<
short_t d,
typename T>
396void gsObjFoldoverFree<d, T>::defaultOptions() {
398 m_options.addReal(
"ff_Delta",
"Sets the delta value", 1e-2);
401template<
short_t d,
typename T>
402void gsObjFoldoverFree<d, T>::addOptions(
const gsOptionList &options) {
403 m_options.update(options, gsOptionList::addIfUnknown);
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);
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);
417 auto EfoldoverFree = (m_eps -
jac(G).det()).ppartval();
418 return m_evaluator.integral(EfoldoverFree);
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);
428 space space1 = m_assembler.getSpace(m_mb, d);
429 space1.setupMapper(m_mapper);
432 auto derJacDet = frprod2(space1,
jac(G).tr().adj());
434 gsConstantFunction<T> zeroFunc(gsVector<T>::Zero(d), d);
435 auto zeroVar = m_evaluator.getVariable(zeroFunc);
437 auto Eder =
ternary(m_eps -
jac(G).det(), -derJacDet, space1 * zeroVar);
439 m_assembler.initSystem();
440 m_assembler.assemble(Eder);
442 result.resize(m_assembler.rhs().rows());
443 std::copy(m_assembler.rhs().data(),
444 m_assembler.rhs().data() + m_assembler.rhs().rows(),
448template<
short_t d,
typename T>
449gsObjQualityImprovePt<d, T>::gsObjQualityImprovePt(
450 const gsMultiPatch<T> &patches,
454 m_mapper(std::move(mapper)),
456 m_assembler.setIntegrationElements(m_mb);
457 m_evaluator = gsExprEvaluator<T>(m_assembler);
461template<
short_t d,
typename T>
462void gsObjQualityImprovePt<d, T>::defaultOptions() {
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);
468template<
short_t d,
typename T>
469void gsObjQualityImprovePt<d, T>::addOptions(
const gsOptionList &options) {
470 m_options.update(options, gsOptionList::addIfUnknown);
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);
481template<
short_t d,
typename T>
482T gsObjQualityImprovePt<d, T>::evalObj(
const gsAsConstVector<T> &u)
const {
483 return evalObj_impl<d>(u);
486template<
short_t d,
typename T>
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);
493 if (m_evaluator.min(
jac(G).det()) < 0) {
494 return std::numeric_limits<T>::max();
496 auto Ewinslow =
jac(G).sqNorm() /
jac(G).det();
497 auto Euniform = pow(
jac(G).det(), 2);
499 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
503template<
short_t d,
typename T>
505typename std::enable_if<_d == 3, T>::type
506gsObjQualityImprovePt<d, T>::evalObj_impl(
const gsAsConstVector<T> &u)
const {
508 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
510 geometryMap G = m_evaluator.getMap(m_mp);
512 if (m_evaluator.min(
jac(G).det()) < 0) {
513 return std::numeric_limits<T>::max();
515 auto Euniform = pow(
jac(G).det(), 2);
516 auto Ewinslow = 0.5 * (
jac(G).sqNorm() *
jac(G).inv().sqNorm());
521 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
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);
531template<
short_t d,
typename T>
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);
538 geometryMap G = m_assembler.getMap(m_mp);
540 space space1 = m_assembler.getSpace(m_mb, d);
541 space1.setupMapper(m_mapper);
544 auto derJacDet = frprod2(space1,
jac(G).tr().adj());
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;
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());
558template<
short_t d,
typename T>
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);
566 space space1 = m_assembler.getSpace(m_mb, d);
567 space1.setupMapper(m_mapper);
570 auto derJacDet = frprod2(space1,
jac(G).tr().adj());
572 auto derEwinslow = frprod2(space1, (
jac(G).inv().sqNorm() *
jac(G) -
jac(G)
573 .sqNorm() * (
jac(G).tr() *
jac(G) *
jac(G).tr()).inv()));
580 auto derEuniform = 2.0 *
jac(G).det() * derJacDet;
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());
589template<
short_t d,
typename T>
590gsObjVHPt<d, T>::gsObjVHPt(
const gsMultiPatch<T> &patches,
594 m_mapper(std::move(mapper)),
597 m_assembler.setIntegrationElements(m_mb);
598 m_evaluator = gsExprEvaluator<T>(m_assembler);
601template<
short_t d,
typename T>
602void gsObjVHPt<d, T>::defaultOptions() {
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);
608template<
short_t d,
typename T>
609void gsObjVHPt<d, T>::addOptions(
const gsOptionList &options) {
610 m_options.update(options, gsOptionList::addIfUnknown);
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);
619template<
short_t d,
typename T>
620T gsObjVHPt<d, T>::evalObj(
const gsAsConstVector<T> &u)
const {
621 return evalObj_impl<d>(u);
624template<
short_t d,
typename T>
626typename std::enable_if<_d == 2, T>::type
627gsObjVHPt<d, T>::evalObj_impl(
const gsAsConstVector<T> &u)
const {
629 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
632 geometryMap G = m_evaluator.getMap(m_mp);
635 space space1 = m_assembler.getSpace(m_mb, d);
636 space1.setupMapper(m_mapper);
639 auto metricMat =
jac(G).tr() *
jac(G);
642 auto hessG = hess(G);
645 auto scale = metricMat.trace();
648 auto Lx = hessG % metricMat.adj() / scale.val();
651 auto nonlinearSystem =
frprod3(space1, Lx);
654 m_assembler.initSystem();
655 m_assembler.assemble(nonlinearSystem);
658 gsVector<T> result = gsAsVector<T>(
659 const_cast<T *
>(m_assembler.rhs().data()),
660 m_assembler.rhs().rows());
662 return result.squaredNorm();
665template<
short_t d,
typename T>
667typename std::enable_if<_d == 3, T>::type
668gsObjVHPt<d, T>::evalObj_impl(
const gsAsConstVector<T> &u)
const {
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);
678template<
short_t d,
typename T>
679gsObjPenaltyPt<d, T>::gsObjPenaltyPt(
const gsMultiPatch<T> &patches,
683 m_mapper(std::move(mapper)),
686 m_assembler.setIntegrationElements(m_mb);
687 m_evaluator = gsExprEvaluator<T>(m_assembler);
690template<
short_t d,
typename T>
691void gsObjPenaltyPt<d, T>::defaultOptions() {
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);
697template<
short_t d,
typename T>
698void gsObjPenaltyPt<d, T>::addOptions(
const gsOptionList &options) {
699 m_options.update(options, gsOptionList::addIfUnknown);
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);
710template<
short_t d,
typename T>
711T gsObjPenaltyPt<d, T>::evalObj(
const gsAsConstVector<T> &u)
const {
713 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
716 geometryMap G = m_evaluator.getMap(m_mp);
719 gsConstantFunction<T> eps1(pow(m_eps, 2.0), d);
720 auto eps = m_evaluator.getVariable(eps1);
724 chi = 0.5 * (
jac(G).det() + pow(eps.val() + pow(
jac(G).det(), 2.0), 0.5));
727 auto Ewinslow =
jac(G).sqNorm() / pow(chi, (2.0 /
static_cast<T
>(d)));
728 auto Euniform = pow(
jac(G).det(), 2.0);
730 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
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);
740template<
short_t d,
typename T>
742typename std::enable_if<_d == 2, T>::type
743gsObjPenaltyPt<d, T>::gradObj_into_impl(
const gsAsConstVector<T> &u,
744 gsAsVector<T> &result)
const {
747 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
750 geometryMap G = m_assembler.getMap(m_mp);
751 space space1 = m_assembler.getSpace(m_mb, d);
752 space1.setupMapper(m_mapper);
755 auto derJacDet =
jac(space1) %
jac(G).tr().adj();
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);
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);
769 auto Ewinslow =
jac(G).sqNorm() / chi;
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;
777 m_assembler.initSystem();
778 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
781 result = gsAsVector<T>(
const_cast<T *
>(m_assembler.rhs().data()),
782 m_assembler.rhs().rows());
787template<
short_t d,
typename T>
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;
795 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
798 geometryMap G = m_assembler.getMap(m_mp);
799 space space1 = m_assembler.getSpace(m_mb, d);
800 space1.setupMapper(m_mapper);
803 auto derJacDet = frprod2(space1,
jac(G).tr().adj());
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);
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);
817 auto Ewinslow =
jac(G).sqNorm() / pow(chi, twoThirds);
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;
825 m_assembler.initSystem();
826 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
829 result = gsAsVector<T>(
const_cast<T *
>(m_assembler.rhs().data()),
830 m_assembler.rhs().rows());
836template<
short_t d,
typename T>
837gsObjPenaltyPt2<d, T>::gsObjPenaltyPt2(
const gsMultiPatch<T> &patches,
841 m_mapper(std::move(mapper)),
844 m_assembler.setIntegrationElements(m_mb);
845 m_evaluator = gsExprEvaluator<T>(m_assembler);
848template<
short_t d,
typename T>
849void gsObjPenaltyPt2<d, T>::defaultOptions() {
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);
855template<
short_t d,
typename T>
856void gsObjPenaltyPt2<d, T>::addOptions(
const gsOptionList &options) {
857 m_options.update(options, gsOptionList::addIfUnknown);
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);
868template<
short_t d,
typename T>
869T gsObjPenaltyPt2<d, T>::evalObj(
const gsAsConstVector<T> &u)
const {
870 return evalObj_impl<d>(u);
873template<
short_t d,
typename T>
875typename std::enable_if<_d == 2, T>::type
876gsObjPenaltyPt2<d, T>::evalObj_impl(
const gsAsConstVector<T> &u)
const {
878 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
881 geometryMap G = m_evaluator.getMap(m_mp);
884 gsConstantFunction<T> epsilon(m_eps, d);
885 auto eps = m_evaluator.getVariable(epsilon);
888 auto chiPPart = eps * ((
jac(G).det() - eps.val()).exp());
890 ternary(eps.val() -
jac(G).det(), chiPPart.val(),
jac(G).det().val());
893 auto Ewinslow =
jac(G).sqNorm() / chi;
895 auto Euniform = pow(
jac(G).det(), 2.0);
898 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
901template<
short_t d,
typename T>
903typename std::enable_if<_d == 3, T>::type
904gsObjPenaltyPt2<d, T>::evalObj_impl(
const gsAsConstVector<T> &u)
const {
906 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
909 geometryMap G = m_evaluator.getMap(m_mp);
912 gsConstantFunction<T> epsilon(m_eps, d);
913 auto eps = m_evaluator.getVariable(epsilon);
916 auto chiPPart = eps * ((
jac(G).det() - eps.val()).exp());
918 ternary(eps.val() -
jac(G).det(), chiPPart.val(),
jac(G).det().val());
921 auto Ewinslow = 0.5 * (
jac(G).sqNorm() *
jac(G).inv().sqNorm()) *
922 pow(
jac(G).det(), 2.0) / pow(chi, 2.0);
926 auto Euniform = pow(
jac(G).det(), 2.0);
929 return m_evaluator.integral(m_lambda1 * Ewinslow + m_lambda2 * Euniform);
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);
938template<
short_t d,
typename T>
940typename std::enable_if<_d == 2, T>::type
941gsObjPenaltyPt2<d, T>::gradObj_into_impl(
const gsAsConstVector<T> &u,
942 gsAsVector<T> &result)
const {
944 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
947 geometryMap G = m_assembler.getMap(m_mp);
950 space spaceMapper = m_assembler.getSpace(m_mb, d);
951 spaceMapper.setupMapper(m_mapper);
954 auto derJacDet = frprod2(spaceMapper,
jac(G).tr().adj());
957 gsConstantFunction<T> epsilon(m_eps, d);
958 gsConstantFunction<T> unity(1.0, d);
961 auto eps = m_evaluator.getVariable(epsilon);
962 auto unit = m_evaluator.getVariable(unity);
965 auto chiPPart = eps * ((
jac(G).det() - eps.val()).exp());
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());
971 auto Ewinslow =
jac(G).sqNorm() / chi;
972 auto derEwinslow = 1.0 / chi * (2.0 * frprod2(spaceMapper,
jac(G)) -
973 Ewinslow.val() * chip * derJacDet);
977 auto derEuniform = 2.0 *
jac(G).det() * derJacDet;
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());
989template<
short_t d,
typename T>
991typename std::enable_if<_d == 3, T>::type
992gsObjPenaltyPt2<d, T>::gradObj_into_impl(
const gsAsConstVector<T> &u,
993 gsAsVector<T> &result)
const {
995 convertFreeVectorToMultiPatch<T>(u, m_mapper, m_mp);
998 geometryMap G = m_assembler.getMap(m_mp);
1001 space space1 = m_assembler.getSpace(m_mb, d);
1002 space1.setupMapper(m_mapper);
1005 auto derJacDet = frprod2(space1,
jac(G).tr().adj());
1008 gsConstantFunction<T> epsilonFunction(m_eps, d);
1009 gsConstantFunction<T> unityFunction(1.0, d);
1012 auto eps = m_evaluator.getVariable(epsilonFunction);
1013 auto unit = m_evaluator.getVariable(unityFunction);
1016 auto chiPPart = eps * ((
jac(G).det() - eps.val()).exp());
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());
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()));
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;
1033 auto derEuniform = 2.0 *
jac(G).det() * derJacDet;
1036 m_assembler.initSystem();
1037 m_assembler.assemble(m_lambda1 * derEwinslow + m_lambda2 * derEuniform);
1040 result = gsAsVector<T>(
const_cast<T *
>(m_assembler.rhs().data()),
1041 m_assembler.rhs().rows());
1042 return EXIT_SUCCESS;
1046template<
short_t d,
typename T>
1066 verboseLog(
"PDE-based parameterization construction ...\n",
1067 options.
askInt(
"Verbose", 0));
1073 assembler.setIntegrationElements(mb);
1077 for (
auto bit = mp.bBegin(); bit != mp.bEnd(); ++bit) {
1081 space space1 = assembler.getSpace(mb, d);
1082 space1.setup(bc, dirichlet::homogeneous, 0);
1085 typedef std::function<gsSparseMatrix<T>(
gsVector<T> const &)> Jacobian_t;
1086 typedef std::function<gsVector<T>(
gsVector<T> const &)> Residual_t;
1089 Residual_t Residual = [&assembler, &mapper, &mpSubstitute, &space1](
1092 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1095 geometryMap G = assembler.getMap(mpSubstitute);
1098 auto metricMat = jac(G).tr() * jac(G);
1099 auto hessG = hess(G);
1102 auto scale = metricMat.trace();
1105 auto Lx = hessG % metricMat.adj() / scale.val();
1108 auto nonlinearSystem = frprod3(space1, Lx);
1111 assembler.initSystem();
1112 assembler.assemble(nonlinearSystem);
1115 return assembler.rhs();
1118 Jacobian_t Jacobian;
1119 int preconditionerType = options.
askInt(
"AAPreconditionType", 0);
1120 switch (preconditionerType) {
1122 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](
gsVector<T> const &x)
1125 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1126 geometryMap G = assembler.getMap(mpSubstitute);
1128 auto jacMat = jacScaledLxDiag(space1, G);
1130 assembler.initSystem();
1131 assembler.assemble(jacMat);
1132 return assembler.matrix();
1137 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](
gsVector<T> const &x)
1140 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1141 geometryMap G = assembler.getMap(mpSubstitute);
1143 auto jacMat = jacScaledLxDiagBlock(space1, G);
1145 assembler.initSystem();
1146 assembler.assemble(jacMat);
1147 return assembler.matrix();
1153 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](
gsVector<T> const &x)
1155 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1156 geometryMap G = assembler.getMap(mpSubstitute);
1158 auto jacMat = jacScaledLx(space1, G);
1160 assembler.initSystem();
1161 assembler.assemble(jacMat);
1162 return assembler.matrix();
1167 param.
m = options.
askInt(
"AAwindowsize", 5);
1175 if (options.
askSwitch(
"needPDEH1",
true))
1177 verboseLog(
"\nStart parameterization improvement by H1 discrezation...",options.
askInt(
"Verbose", 0));
1178 Residual = [&assembler, &mapper, &mpSubstitute, &space1](
gsVector<T> const &x)
1181 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1182 geometryMap G = assembler.getMap(mpSubstitute);
1184 auto invJacMat = jac(G).inv();
1185 auto jacU = jac(space1) * invJacMat;
1186 auto nonlinearSystem = jacU % invJacMat;
1188 assembler.initSystem();
1189 assembler.assemble(nonlinearSystem);
1190 return assembler.rhs();
1193 preconditionerType = options.
askInt(
"AAPreconditionType", 0);
1194 switch (preconditionerType) {
1196 Jacobian = [&assembler, &mapper, &mpSubstitute, &space1](
gsVector<T> const &x)
1199 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1200 geometryMap G = assembler.getMap(mpSubstitute);
1202 auto jacMat = jacScaledLxH1DiagBlock(space1, G);
1204 assembler.initSystem();
1205 assembler.assemble(jacMat);
1207 return assembler.matrix();
1212 Jacobian = [&assembler, &mapper, &mpSubstitute,&space1](
gsVector<T> const &x)
1215 convertFreeVectorToMultiPatch<T>(x, mapper, mpSubstitute);
1216 geometryMap G = assembler.getMap(mpSubstitute);
1218 auto jacMat = jacScaledLxH1(space1, G);
1220 assembler.initSystem();
1221 assembler.assemble(jacMat);
1223 return assembler.matrix();
1228 solVector = AASolver2.
compute(solVector, Residual, Jacobian);
1232 convertFreeVectorToMultiPatch<T>(solVector, mapper, result);
1241template<
typename E1,
typename E2>
1243template<
typename E1,
typename E2>
1245template<
class E0,
class E1,
class E2>
1248class jacScaledLx_expr;
1250class jacScaledLxDiag_expr;
1252class jacScaledLxDiagBlock_expr;
1254class jacScaledLxH1_expr;
1256class jacScaledLxH1DiagBlock_expr;
1283template<
typename E1,
typename E2>
1284class frprod3_expr :
public _expr<frprod3_expr<E1, E2> > {
1286 typedef typename E2::Scalar Scalar;
1287 enum { ScalarValued = 0, Space = E1::Space, ColBlocks = 0 };
1290 typename E1::Nested_t _u;
1291 typename E2::Nested_t _v;
1293 mutable gsMatrix<Scalar> res, bGrads, b;
1297 frprod3_expr(_expr<E1>
const &u, _expr<E2>
const &v)
1306 const gsMatrix<Scalar> &eval(
const index_t k)
const
1308 auto A = _v.eval(k);
1312 res.noalias() = b * A;
1316 index_t rows()
const {
return 1; }
1317 index_t cols()
const {
return 1; }
1319 void parse(gsExprHelper<Scalar> &evList)
const {
1326 const gsFeSpace<Scalar> &rowVar()
const {
return _u.rowVar(); }
1327 const gsFeSpace<Scalar> &colVar()
const {
return _v.rowVar(); }
1329 void print(std::ostream &os)
const {
1339template<
typename E1,
typename E2>
1343 return frprod3_expr<E1, E2>(u, M);
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;
1352 typedef typename E1::Scalar Scalar;
1354 explicit ternary_expr(_expr<E0>
const &u,
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).");
1367 "Both v and w must be space (or not), but E1::Space = "
1368 << E1::Space <<
" and E2::Space = " << E2::Space);
1370 "Rows of v and w differ. _v.rows() = " << _v.rows()
1374 "Columns of v and w differ. _v.cols() = " << _v.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.");
1382 ScalarValued = E1::ScalarValued,
1383 ColBlocks = E1::ColBlocks,
1390 const Temporary_t eval(
const index_t k)
const {
1391 return (_u.eval(k) > 0 ? _v.eval(k) : _w.eval(k));
1396 index_t rows()
const {
return _v.rows(); }
1397 index_t cols()
const {
return _v.cols(); }
1398 void parse(gsExprHelper<Scalar> &evList)
const {
1404 const gsFeSpace<Scalar> &rowVar()
const {
return _v.rowVar(); }
1405 const gsFeSpace<Scalar> &colVar()
const {
return _v.colVar(); }
1407 void print(std::ostream &os)
const { os <<
"ternary"; }
1453class jacScaledLx_expr :
public _expr<jacScaledLx_expr<E> > {
1455 typedef typename E::Scalar Scalar;
1458 typename E::Nested_t _u;
1459 typename gsGeometryMap<Scalar>::Nested_t _G;
1462 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1464 jacScaledLx_expr(
const E &u,
const gsGeometryMap<Scalar> &G) : _u(u), _G(G) {}
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;
1472 const gsMatrix<Scalar> &eval(
const index_t k)
const {
1473 gsMatrix<Scalar> basis = _u.data().values[0].col(k);
1475 derivBasis = _u.data().values[1].col(k).transpose();
1476 deriv2Basis = _u.data().values[2].col(k).transpose();
1479 deriv2Basis.blockTransposeInPlace(1 + _u.dim());
1481 derivGeom = _G.data().values[1].col(k);
1482 deriv2Geom = _G.data().values[2].col(k);
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);
1488 Scalar scaleFactor = g11 + g22;
1491 (g22 * deriv2Geom(0) + g11 * deriv2Geom(1) - 2.0 * g12 * deriv2Geom(2))
1494 (g22 * deriv2Geom(3) + g11 * deriv2Geom(4) - 2.0 * g12 * deriv2Geom(5))
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);
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;
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;
1522 const index_t A = _u.cardinality() / _u.dim();
1523 res.resize(_u.cardinality(), _u.cardinality());
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;
1533 index_t rows()
const {
return 1; }
1535 index_t cols()
const {
return 1; }
1537 void parse(gsExprHelper<Scalar> &evList)
const {
1545 const gsFeSpace<Scalar> &rowVar()
const {
return _u.rowVar(); }
1546 const gsFeSpace<Scalar> &colVar()
const {
return _u.rowVar(); }
1548 index_t cardinality_impl()
const {
return _u.cardinality_impl(); }
1550 void print(std::ostream &os)
const {
1551 os <<
"jacScaledLx(";
1561class jacScaledLxDiag_expr :
public _expr<jacScaledLxDiag_expr<E> > {
1563 typedef typename E::Scalar Scalar;
1566 typename E::Nested_t _u;
1567 typename gsGeometryMap<Scalar>::Nested_t _G;
1570 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1572 jacScaledLxDiag_expr(
const E &u,
const gsGeometryMap<Scalar> &G)
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;
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();
1586 deriv2Basis.blockTransposeInPlace(1 + _u.dim());
1588 derivGeom = _G.data().values[1].col(k);
1589 deriv2Geom = _G.data().values[2].col(k);
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);
1595 Scalar scaleFactor = g11 + g22;
1598 (g22 * deriv2Geom(0) + g11 * deriv2Geom(1) - 2.0 * g12 * deriv2Geom(2))
1601 (g22 * deriv2Geom(3) + g11 * deriv2Geom(4) - 2.0 * g12 * deriv2Geom(5))
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);
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;
1620 dLxdx = (dLxdx - Lx * (dg11dx + dg22dx)) / scaleFactor;
1621 dLydy = (dLydy - Ly * (dg11dy + dg22dy)) / scaleFactor;
1623 const index_t A = _u.cardinality() / _u.dim();
1624 res.resize(_u.cardinality(), _u.cardinality());
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();
1634 index_t rows()
const {
return 1; }
1636 index_t cols()
const {
return 1; }
1638 void parse(gsExprHelper<Scalar> &evList)
const {
1646 const gsFeSpace<Scalar> &rowVar()
const {
return _u.rowVar(); }
1647 const gsFeSpace<Scalar> &colVar()
const {
return _u.rowVar(); }
1649 index_t cardinality_impl()
const {
return _u.cardinality_impl(); }
1651 void print(std::ostream &os)
const {
1652 os <<
"jacScaledLxDiag(";
1662class jacScaledLxDiagBlock_expr :
public _expr<jacScaledLxDiagBlock_expr<E> > {
1664 typedef typename E::Scalar Scalar;
1667 typename E::Nested_t _u;
1668 typename gsGeometryMap<Scalar>::Nested_t _G;
1671 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1673 jacScaledLxDiagBlock_expr(
const E &u,
const gsGeometryMap<Scalar> &G)
1676 mutable gsMatrix<Scalar> res, basis, derivGeom, deriv2Geom,
derivBasis,
1678 mutable gsMatrix<Scalar> dg11dx, dg11dy, dg22dx, dg22dy, dg12dx, dg12dy;
1679 mutable gsMatrix<Scalar> commonTerm, dLxdx, dLydy;
1683 const gsMatrix<Scalar> &eval(
const index_t k)
const {
1686 derivBasis = _u.data().values[1].col(k).transpose();
1687 deriv2Basis = _u.data().values[2].col(k).transpose();
1690 deriv2Basis.blockTransposeInPlace(1 + _u.dim());
1692 derivGeom = _G.data().values[1].col(k);
1693 deriv2Geom = _G.data().values[2].col(k);
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;
1701 (g22 * deriv2Geom(0) + g11 * deriv2Geom(1) - 2.0 * g12 * deriv2Geom(2))
1704 (g22 * deriv2Geom(3) + g11 * deriv2Geom(4) - 2.0 * g12 * deriv2Geom(5))
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);
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;
1723 dLxdx = (dLxdx - Lx * (dg11dx + dg22dx)) / scaleFactor;
1724 dLydy = (dLydy - Ly * (dg11dy + dg22dy)) / scaleFactor;
1726 const index_t A = _u.cardinality() / _u.dim();
1727 res.resize(_u.cardinality(), _u.cardinality());
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;
1737 index_t rows()
const {
return 1; }
1739 index_t cols()
const {
return 1; }
1741 void parse(gsExprHelper<Scalar> &evList)
const {
1749 const gsFeSpace<Scalar> &rowVar()
const {
return _u.rowVar(); }
1750 const gsFeSpace<Scalar> &colVar()
const {
return _u.rowVar(); }
1752 index_t cardinality_impl()
const {
return _u.cardinality_impl(); }
1754 void print(std::ostream &os)
const {
1755 os <<
"jacScaledLxDiagBlock(";
1765class jacScaledLxH1_expr :
public _expr<jacScaledLxH1_expr<E> > {
1767 typedef typename E::Scalar Scalar;
1770 typename E::Nested_t _u;
1771 typename gsGeometryMap<Scalar>::Nested_t _G;
1774 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1776 jacScaledLxH1_expr(
const E &u,
const gsGeometryMap<Scalar> &G)
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();
1789 derivBasis = _u.data().values[1].col(k).transpose();
1791 gsMatrix<Scalar> invJacHatG;
1792 invJacHatG.noalias() = invJacMat *
derivBasis;
1794 const index_t N = _u.cardinality() / _u.dim();
1795 gsMatrix<Scalar> jacdLxdx(N, N);
1796 gsMatrix<Scalar> jacdLxdy(N, N);
1797 gsMatrix<Scalar> jacdLydx(N, N);
1798 gsMatrix<Scalar> jacdLydy(N, N);
1800 gsMatrix<> temp(2, N);
1801 for (
auto i = 0; i < N; ++i) {
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);
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);
1825 res.resize(_u.cardinality(), _u.cardinality());
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;
1835 index_t rows()
const {
return 1; }
1837 index_t cols()
const {
return 1; }
1839 void parse(gsExprHelper<Scalar> &evList)
const {
1841 _u.data().flags |= NEED_GRAD;
1847 const gsFeSpace<Scalar> &rowVar()
const {
return _u.rowVar(); }
1848 const gsFeSpace<Scalar> &colVar()
const {
return _u.rowVar(); }
1850 index_t cardinality_impl()
const {
return _u.cardinality_impl(); }
1852 void print(std::ostream &os)
const {
1853 os <<
"jacScaledLx(";
1863class jacScaledLxH1DiagBlock_expr
1864 :
public _expr<jacScaledLxH1DiagBlock_expr<E> > {
1866 typedef typename E::Scalar Scalar;
1869 typename E::Nested_t _u;
1870 typename gsGeometryMap<Scalar>::Nested_t _G;
1873 enum { Space = 3, ScalarValued = 0, ColBlocks = 0 };
1875 jacScaledLxH1DiagBlock_expr(
const E &u,
const gsGeometryMap<Scalar> &G)
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();
1888 derivBasis = _u.data().values[1].col(k).transpose();
1890 gsMatrix<Scalar> invJacHatG = invJacMat *
derivBasis;
1892 const index_t N = _u.cardinality() / _u.dim();
1893 gsMatrix<Scalar> jacdLxdx(N, N);
1896 gsMatrix<Scalar> jacdLydy(N, N);
1898 gsMatrix<> temp(2, N);
1899 for (
auto i = 0; i < N; ++i) {
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);
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);
1917 jacdLydy.col(i) = -temp.transpose() * invJacMat.col(1)
1918 - invJacHatG.transpose() * dinvJacdy.col(1);
1921 res.resize(_u.cardinality(), _u.cardinality());
1923 res.topLeftCorner(N, N) = jacdLxdx;
1926 res.bottomRightCorner(N, N) = jacdLydy;
1931 index_t rows()
const {
return 1; }
1933 index_t cols()
const {
return 1; }
1935 void parse(gsExprHelper<Scalar> &evList)
const {
1937 _u.data().flags |= NEED_GRAD;
1943 const gsFeSpace<Scalar> &rowVar()
const {
return _u.rowVar(); }
1944 const gsFeSpace<Scalar> &colVar()
const {
return _u.rowVar(); }
1946 index_t cardinality_impl()
const {
return _u.cardinality_impl(); }
1948 void print(std::ostream &os)
const {
1949 os <<
"jacScaledLx(";
1980template<
typename E1,
typename E2>
1981class frprod2_expr :
public _expr<frprod2_expr<E1, E2> > {
1983 typedef typename E2::Scalar Scalar;
1984 enum { ScalarValued = 0, Space = E1::Space, ColBlocks = 0 };
1987 typename E1::Nested_t _u;
1988 typename E2::Nested_t _v;
1990 mutable gsMatrix<Scalar> res, bGrads;
1994 frprod2_expr(_expr<E1>
const &u, _expr<E2>
const &v)
2003 const gsMatrix<Scalar> &eval(
const index_t k)
const
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());
2012 res.transposeInPlace();
2016 index_t rows()
const {
return 1; }
2017 index_t cols()
const {
return 1; }
2019 void parse(gsExprHelper<Scalar> &evList)
const {
2022 _u.data().flags |= NEED_GRAD;
2025 const gsFeSpace<Scalar> &rowVar()
const {
return _u.rowVar(); }
2026 const gsFeSpace<Scalar> &colVar()
const {
return _v.rowVar(); }
2028 void print(std::ostream &os)
const {
2041 const gsGeometryMap<typename E::Scalar> &G) {
2042 return jacScaledLx_expr<E>(u, G);
2049 const gsGeometryMap<typename E::Scalar> &G) {
2050 return jacScaledLxDiag_expr<E>(u, G);
2057 const gsGeometryMap<typename E::Scalar> &G) {
2058 return jacScaledLxDiagBlock_expr<E>(u, G);
2065 const gsGeometryMap<typename E::Scalar> &G) {
2066 return jacScaledLxH1_expr<E>(u, G);
2073 const gsGeometryMap<
2074 typename E::Scalar> &G) {
2075 return jacScaledLxH1DiagBlock_expr<E>(u, G);
2079template<
typename E1,
typename E2>
2081frprod2_expr<E1, E2>
const frprod2(E1
const &u,
2083 return frprod2_expr<E1, E2>(u, M);
2087template<
class E0,
class E1,
class E2>
2092 return ternary_expr<E0, E1, E2>(u, v, w);
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