G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsBiharmonicExprAssembler.hpp
1
16#pragma once
17
23#include <gsCore/gsBasis.h>
25
26#ifdef gsSpectra_ENABLED
27#include <gsSpectra/gsSpectra.h>
28#endif
29
30namespace gismo
31{
32
33template <class T>
34gsBiharmonicExprAssembler<T>::gsBiharmonicExprAssembler(const gsMultiPatch<T> & mp,
35 const gsMultiBasis<T> & mb,
36 const gsFunctionSet<T> & force,
37 const gsBoundaryConditions<T> & bcs
38 )
39:
40m_patches(mp),
41m_basis(mb),
42m_spaceBasis(&mb),
43m_force(&force),
44m_bcs(bcs)
45{
46 this->_defaultOptions();
47 this->_getOptions();
48 this->_initialize();
49};
50
51template <class T>
52gsBiharmonicExprAssembler<T>& gsBiharmonicExprAssembler<T>::operator=( const gsBiharmonicExprAssembler& other )
53{
54 if (this!=&other)
55 {
56 m_penalty=other.m_penalty;
57 m_lambda=other.m_lambda;
58 m_second=other.m_second;
59 m_continuity=other.m_continuity;
60
61 m_patches=other.m_patches;
62 m_basis=other.m_basis;
63 m_spaceBasis=other.m_spaceBasis;
64 m_bcs=other.m_bcs;
65 m_force=other.m_force;
66
67 m_options=other.m_options;
68
69 // To do: make copy constructor for the gsExprAssembler
70 m_assembler.setIntegrationElements(m_basis);
71 m_assembler.setOptions(m_options);
72 }
73 return *this;
74}
75
76template <class T>
77gsBiharmonicExprAssembler<T>& gsBiharmonicExprAssembler<T>::operator=( gsBiharmonicExprAssembler&& other )
78{
79 m_penalty=give(other.m_penalty);
80 m_lambda=give(other.m_lambda);
81 m_second=give(other.m_second);
82 m_continuity=give(other.m_continuity);
83
84 m_patches=give(other.m_patches);
85 m_basis=give(other.m_basis);
86 m_spaceBasis=give(other.m_spaceBasis);
87 m_bcs=give(other.m_bcs);
88 m_force=give(other.m_force);
89
90 m_options=give(other.m_options);
91
92 // To do: make copy constructor for the gsExprAssembler
93 m_assembler.setIntegrationElements(m_basis);
94 m_assembler.setOptions(m_options);
95 return *this;
96}
97
98template <class T>
99void gsBiharmonicExprAssembler<T>::_defaultOptions()
100{
101 m_options.addReal("PenaltyIfc","Parameter Nitsche's method",-1);
102 m_options.addReal("Lambda","Parameter for BC projection",1e-5);
103 m_options.addSwitch("Second","Assemble the second biharmonic equation",false);
104 m_options.addInt("Continuity","Set the continuity for the space",-1);
105 // Assembler options
106 gsOptionList assemblerOptions = m_assembler.defaultOptions().wrapIntoGroup("ExprAssembler");
107 m_options.update(assemblerOptions,gsOptionList::addIfUnknown);
108}
109
110template <class T>
111void gsBiharmonicExprAssembler<T>::_getOptions()
112{
113 m_penalty = m_options.getReal("PenaltyIfc");
114 m_lambda = m_options.getReal("Lambda");
115 m_second = m_options.getSwitch("Second");
116 m_continuity = m_options.getInt("Continuity");
117
118 GISMO_ENSURE(m_options.hasGroup("ExprAssembler"),"The option list does not contain options with the label 'ExprAssembler'!");
119 m_assembler.setOptions(m_options.getGroup("ExprAssembler"));
120}
121
122template<class T>
123void gsBiharmonicExprAssembler<T>::setOptions(gsOptionList & options)
124{
125 m_options.update(options,gsOptionList::ignoreIfUnknown);
126 this->_getOptions();
127 this->_initialize();
128}
129
130template <class T>
131void gsBiharmonicExprAssembler<T>::_initialize()
132{
133 // Elements used for numerical integration
134 m_assembler.setIntegrationElements(m_basis);
135 m_assembler.setOptions(m_options.getGroup("ExprAssembler"));
136
137 GISMO_ASSERT(m_bcs.hasGeoMap(),"No geometry map was assigned to the boundary conditions. Use bc.setGeoMap to assign one!");
138
139 m_assembler.getSpace(*m_spaceBasis, 1, 0); // last argument is the space ID
140}
141
142template <class T>
143void gsBiharmonicExprAssembler<T>::_setup(const expr::gsFeSpace<T> & u)
144{
145 // Setup the mapper
146 gsDofMapper map;
147 _setMapperForBiharmonic(m_bcs, *m_spaceBasis, map);
148
149 // Setup the system
150 u.setupMapper(map);
151 _getDirichletNeumannValuesL2Projection(m_patches, m_basis, m_bcs, *m_spaceBasis, u);
152}
153
154template <class T>
155void gsBiharmonicExprAssembler<T>::assembleMass()
156{
157 m_assembler.cleanUp();
158 this->_getOptions();
159 this->_initialize();
160
161 // Set the geometry map
162 geometryMap G = m_assembler.getMap(m_patches);
163
164 // Set the discretization space
165 auto u = m_assembler.trialSpace(0);
166
167 // Setup the mapper
168 this->_setup(u);
169
170 // Initialize the system
171 m_assembler.initSystem();
172
173 // Compute the system matrix and right-hand side
174 m_assembler.assemble(u * u.tr() * meas(G));
175}
176
177template <class T>
178void gsBiharmonicExprAssembler<T>::assemble()
179{
180 m_assembler.cleanUp();
181 this->_getOptions();
182 this->_initialize();
183
184 // Set the geometry map
185 geometryMap G = m_assembler.getMap(m_patches);
186
187 // Set the source term
188 auto ff = m_assembler.getCoeff(*m_force, G);
189
190 // Set the discretization space
191 auto u = m_assembler.trialSpace(0);
192
193 // Setup the mapper
194 this->_setup(u);
195
196 // Initialize the system
197 m_assembler.initSystem();
198
199 // Compute the system matrix and right-hand side
200 m_assembler.assemble(ilapl(u, G) * ilapl(u, G).tr() * meas(G),u * ff * meas(G));
201
202 // Enforce Laplace conditions to right-hand side
203 auto g_L = m_assembler.getBdrFunction(G); // Set the laplace bdy value
204 //auto g_L = A.getCoeff(laplace, G);
205 m_assembler.assembleBdr(m_bcs.get("Laplace"), (igrad(u, G) * nv(G)) * g_L.tr() );
206
207 // Optionally, assemble Nitsche
208 if (m_continuity>0)
209 {
210 gsMatrix<T> mu_interfaces;
211 if (m_penalty == -1)
212 _computeStabilityParameter(m_patches, m_basis, mu_interfaces);
213
214 index_t i = 0;
215 for ( typename gsMultiPatch<T>::const_iiterator it = m_patches.iBegin(); it != m_patches.iEnd(); ++it, ++i)
216 {
217 T stab = 4 * ( m_basis.maxCwiseDegree() + m_basis.dim() ) * ( m_basis.maxCwiseDegree() + 1 );
218 T m_h = m_basis.basis(0).getMinCellLength(); // m_basis.basis(0).getMinCellLength();
219 T mu = 2 * stab / m_h;
220 T alpha = 1;
221
222 //mu = penalty_init == -1.0 ? mu : penalty_init / m_h;
223 if (m_penalty == -1)
224 mu = mu_interfaces(i,0) / m_h;
225 else
226 mu = m_penalty / m_h;
227
228 std::vector<boundaryInterface> iFace;
229 iFace.push_back(*it);
230 m_assembler.assembleIfc(iFace,
231 //B11
232 -alpha * 0.5 * igrad(u.left(), G) * nv(G.left()).normalized() *
233 (ilapl(u.left(), G)).tr() * nv(G.left()).norm(),
234 -alpha * 0.5 *
235 (igrad(u.left(), G) * nv(G.left()).normalized() * (ilapl(u.left(), G)).tr()).tr() *
236 nv(G.left()).norm(),
237 //B12
238 -alpha * 0.5 * igrad(u.left(), G.left()) * nv(G.left()).normalized() *
239 (ilapl(u.right(), G.right())).tr() * nv(G.left()).norm(),
240 -alpha * 0.5 * (igrad(u.left(), G.left()) * nv(G.left()).normalized() *
241 (ilapl(u.right(), G.right())).tr()).tr() * nv(G.left()).norm(),
242 //B21
243 alpha * 0.5 * igrad(u.right(), G.right()) * nv(G.left()).normalized() *
244 (ilapl(u.left(), G.left())).tr() * nv(G.left()).norm(),
245 alpha * 0.5 * (igrad(u.right(), G.right()) * nv(G.left()).normalized() *
246 (ilapl(u.left(), G.left())).tr()).tr() * nv(G.left()).norm(),
247 //B22
248 alpha * 0.5 * igrad(u.right(), G.right()) * nv(G.left()).normalized() *
249 (ilapl(u.right(), G.right())).tr() * nv(G.left()).norm(),
250 alpha * 0.5 * (igrad(u.right(), G.right()) * nv(G.left()).normalized() *
251 (ilapl(u.right(), G.right())).tr()).tr() * nv(G.left()).norm(),
252
253 // E11
254 mu * igrad(u.left(), G.left()) * nv(G.left()).normalized() *
255 (igrad(u.left(), G.left()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm(),
256 //-E12
257 -mu * (igrad(u.left(), G.left()) * nv(G.left()).normalized()) *
258 (igrad(u.right(), G.right()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm(),
259 //-E21
260 -mu * (igrad(u.right(), G.right()) * nv(G.left()).normalized()) *
261 (igrad(u.left(), G.left()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm(),
262 // E22
263 mu * igrad(u.right(), G.right()) * nv(G.left()).normalized() *
264 (igrad(u.right(), G.right()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm()
265 );
266
267 if (m_continuity>1)
268 gsDebug<<"Put here the terms for continuity 1\n";
269 }
270 }
271}
272
273
274template <class T>
275void gsBiharmonicExprAssembler<T>::assembleRHS()
276{
277 m_assembler.cleanUp();
278 this->_getOptions();
279 this->_initialize();
280
281 // Set the geometry map
282 geometryMap G = m_assembler.getMap(m_patches);
283
284 // Set the source term
285 auto ff = m_assembler.getCoeff(*m_force, G);
286
287 // Set the discretization space
288 auto u = m_assembler.trialSpace(0);
289
290 // Setup the mapper
291 this->_setup(u);
292
293 // Initialize the system
294 m_assembler.initVector(1);
295
296 // Compute the system matrix and right-hand side
297 m_assembler.assemble(u * ff * meas(G));
298
299 // Enforce Laplace conditions to right-hand side
300 auto g_L = m_assembler.getBdrFunction(G); // Set the laplace bdy value
301 //auto g_L = A.getCoeff(laplace, G);
302 m_assembler.assembleBdr(m_bcs.get("Laplace"), (igrad(u, G) * nv(G)) * g_L.tr() );
303}
304
305template <class T>
306void gsBiharmonicExprAssembler<T>::assembleLHS()
307{
308 m_assembler.cleanUp();
309 this->_getOptions();
310 this->_initialize();
311
312 // Set the geometry map
313 geometryMap G = m_assembler.getMap(m_patches);
314
315 // Set the discretization space
316 auto u = m_assembler.trialSpace(0);
317
318 // Setup the mapper
319 this->_setup(u);
320
321 // Initialize the system
322 m_assembler.initMatrix();
323
324 // Compute the system matrix and right-hand side
325 m_assembler.assemble(ilapl(u, G) * ilapl(u, G).tr() * meas(G));
326
327 // Optionally, assemble Nitsche
328 if (m_continuity>0)
329 {
330 gsMatrix<T> mu_interfaces;
331 if (m_penalty == -1)
332 _computeStabilityParameter(m_patches, m_basis, mu_interfaces);
333
334 index_t i = 0;
335 for ( typename gsMultiPatch<T>::const_iiterator it = m_patches.iBegin(); it != m_patches.iEnd(); ++it, ++i)
336 {
337 T stab = 4 * ( m_basis.maxCwiseDegree() + m_basis.dim() ) * ( m_basis.maxCwiseDegree() + 1 );
338 T m_h = m_basis.basis(0).getMinCellLength(); // m_basis.basis(0).getMinCellLength();
339 T mu = 2 * stab / m_h;
340 T alpha = 1;
341
342 //mu = penalty_init == -1.0 ? mu : penalty_init / m_h;
343 if (m_penalty == -1)
344 mu = mu_interfaces(i,0) / m_h;
345 else
346 mu = m_penalty / m_h;
347
348 std::vector<boundaryInterface> iFace;
349 iFace.push_back(*it);
350 m_assembler.assembleIfc(iFace,
351 //B11
352 -alpha * 0.5 * igrad(u.left(), G) * nv(G.left()).normalized() *
353 (ilapl(u.left(), G)).tr() * nv(G.left()).norm(),
354 -alpha * 0.5 *
355 (igrad(u.left(), G) * nv(G.left()).normalized() * (ilapl(u.left(), G)).tr()).tr() *
356 nv(G.left()).norm(),
357 //B12
358 -alpha * 0.5 * igrad(u.left(), G.left()) * nv(G.left()).normalized() *
359 (ilapl(u.right(), G.right())).tr() * nv(G.left()).norm(),
360 -alpha * 0.5 * (igrad(u.left(), G.left()) * nv(G.left()).normalized() *
361 (ilapl(u.right(), G.right())).tr()).tr() * nv(G.left()).norm(),
362 //B21
363 alpha * 0.5 * igrad(u.right(), G.right()) * nv(G.left()).normalized() *
364 (ilapl(u.left(), G.left())).tr() * nv(G.left()).norm(),
365 alpha * 0.5 * (igrad(u.right(), G.right()) * nv(G.left()).normalized() *
366 (ilapl(u.left(), G.left())).tr()).tr() * nv(G.left()).norm(),
367 //B22
368 alpha * 0.5 * igrad(u.right(), G.right()) * nv(G.left()).normalized() *
369 (ilapl(u.right(), G.right())).tr() * nv(G.left()).norm(),
370 alpha * 0.5 * (igrad(u.right(), G.right()) * nv(G.left()).normalized() *
371 (ilapl(u.right(), G.right())).tr()).tr() * nv(G.left()).norm(),
372
373 // E11
374 mu * igrad(u.left(), G.left()) * nv(G.left()).normalized() *
375 (igrad(u.left(), G.left()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm(),
376 //-E12
377 -mu * (igrad(u.left(), G.left()) * nv(G.left()).normalized()) *
378 (igrad(u.right(), G.right()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm(),
379 //-E21
380 -mu * (igrad(u.right(), G.right()) * nv(G.left()).normalized()) *
381 (igrad(u.left(), G.left()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm(),
382 // E22
383 mu * igrad(u.right(), G.right()) * nv(G.left()).normalized() *
384 (igrad(u.right(), G.right()) * nv(G.left()).normalized()).tr() * nv(G.left()).norm()
385 );
386
387 if (m_continuity>1)
388 gsDebug<<"Put here the terms for continuity 1\n";
389 }
390 }
391}
392
393template <class T>
394void gsBiharmonicExprAssembler<T>::constructSolution(gsMatrix<T> & solVector)
395{
396 auto u = m_assembler.trialSpace(0);
397 auto u_sol = m_assembler.getSolution(u, solVector);
398
399 if (const gsMappedBasis<2,T> * bb2 = dynamic_cast<const gsMappedBasis<2,T> *>(m_spaceBasis))
400 {
401 u_sol.extract(m_mspline);
402 }
403 else
404 {
405 u_sol.extract(m_sol);
406 }
407
408}
409
410template <class T>
411typename gsFunctionSet<T>::Ptr gsBiharmonicExprAssembler<T>::getSolution() const
412{
413 if (const gsMappedBasis<2,T> * bb2 = dynamic_cast<const gsMappedBasis<2,T> *>(m_spaceBasis))
414 {
415 return memory::make_shared_not_owned(static_cast<gsFunctionSet<T>*>(&m_mspline));
416 }
417 else
418 {
419 return memory::make_shared_not_owned(static_cast<gsFunctionSet<T>*>(&m_sol));
420 }
421}
422
423template <class T>
424T gsBiharmonicExprAssembler<T>::l2error(gsMatrix<T> & solVector, const gsFunctionSet<T> & exact)
425{
426 auto u = m_assembler.trialSpace(0);
427 auto u_sol = m_assembler.getSolution(u, solVector);
428 geometryMap G = m_assembler.getMap(m_patches);
429
430 gsExprEvaluator<T> ev(m_assembler);
431 auto u_ex = ev.getVariable(exact, G);
432
433 T l2err= math::sqrt( ev.integral( (u_ex - u_sol).sqNorm() * meas(G) ) );
434 return l2err;
435}
436
437template <class T>
438T gsBiharmonicExprAssembler<T>::h1error(gsMatrix<T> & solVector, const gsFunctionSet<T> & exact)
439{
440 geometryMap G = m_assembler.getMap(m_patches);
441 auto u = m_assembler.trialSpace(0);
442 auto u_sol = m_assembler.getSolution(u, solVector);
443
444 gsExprEvaluator<T> ev(m_assembler);
445 auto u_ex = ev.getVariable(exact, G);
446
447 T l2err = math::sqrt( ev.integral( (u_ex - u_sol).sqNorm() * meas(G) ) ); // / ev.integral(f.sqNorm()*meas(G)) );
448 T h1err = l2err +
449 math::sqrt(ev.integral( ( igrad(u_ex) - igrad(u_sol,G) ).sqNorm() * meas(G) ));
450 return h1err;
451}
452
453template <class T>
454T gsBiharmonicExprAssembler<T>::h2error(gsMatrix<T> & solVector, const gsFunctionSet<T> & exact)
455{
456 geometryMap G = m_assembler.getMap(m_patches);
457 auto u = m_assembler.trialSpace(0);
458 auto u_sol = m_assembler.getSolution(u, solVector);
459
460 gsExprEvaluator<T> ev(m_assembler);
461 auto u_ex = ev.getVariable(exact, G);
462
463 T l2err = math::sqrt( ev.integral( (u_ex - u_sol).sqNorm() * meas(G) ) ); // / ev.integral(f.sqNorm()*meas(G)) );
464 T h1err = l2err +
465 math::sqrt(ev.integral( ( igrad(u_ex) - igrad(u_sol,G) ).sqNorm() * meas(G) ));
466 T h2err = h1err +
467 math::sqrt(ev.integral( ( ihess(u_ex) - ihess(u_sol,G) ).sqNorm() * meas(G) )); // /ev.integral( ihess(f).sqNorm()*meas(G) )
468 return h2err;
469}
470
471
472template <class T>
473std::tuple<T,T,T> gsBiharmonicExprAssembler<T>::errors(gsMatrix<T> & solVector, const gsFunctionSet<T> & exact)
474{
475 geometryMap G = m_assembler.getMap(m_patches);
476 auto u = m_assembler.trialSpace(0);
477 auto u_sol = m_assembler.getSolution(u, solVector);
478
479 gsExprEvaluator<T> ev(m_assembler);
480 auto u_ex = ev.getVariable(exact, G);
481
482 T l2err = math::sqrt( ev.integral( (u_ex - u_sol).sqNorm() * meas(G) ) ); // / ev.integral(f.sqNorm()*meas(G)) );
483 T h1err = l2err +
484 math::sqrt(ev.integral( ( igrad(u_ex) - igrad(u_sol,G) ).sqNorm() * meas(G) ));
485 T h2err = h1err +
486 math::sqrt(ev.integral( ( ihess(u_ex) - ihess(u_sol,G) ).sqNorm() * meas(G) )); // /ev.integral( ihess(f).sqNorm()*meas(G) )
487 return std::make_tuple(l2err,h1err,h2err);
488}
489
490template <class T>
491T gsBiharmonicExprAssembler<T>::interfaceError(gsMatrix<T> & solVector, const gsFunctionSet<T> & exact)
492{
493 GISMO_UNUSED(exact);
494 if (const gsMappedBasis<2,T> * bb2 = dynamic_cast<const gsMappedBasis<2,T> *>(m_spaceBasis))
495 {
496 geometryMap G = m_assembler.getMap(m_patches);
497 auto u = m_assembler.trialSpace(0);
498 auto u_sol = m_assembler.getSolution(u, solVector);
499
500 gsMatrix<T> solFull;
501 u_sol.extractFull(solFull);
502 gsMappedSpline<2, T> mappedSpline(*bb2, solFull);
503
504 auto ms_sol = m_assembler.getCoeff(mappedSpline);
505
506 gsExprEvaluator<T> ev(m_assembler);
507 T IFaceErr = math::sqrt(ev.integralInterface(((igrad(ms_sol.left(), G.left()) -
508 igrad(ms_sol.right(), G.right())) *
509 nv(G).normalized()).sqNorm() * meas(G)));
510
511 return IFaceErr;
512 }
513 else
514 {
515 geometryMap G = m_assembler.getMap(m_patches);
516 auto u = m_assembler.trialSpace(0);
517 auto u_sol = m_assembler.getSolution(u, solVector);
518
519 gsMultiPatch<T> sol;
520 u_sol.extract(sol);
521
522 auto ms_sol = m_assembler.getCoeff(sol);
523
524 gsExprEvaluator<T> ev(m_assembler);
525 T IFaceErr = math::sqrt(ev.integralInterface(((igrad(ms_sol.left(), G.left()) -
526 igrad(ms_sol.right(), G.right())) *
527 nv(G).normalized()).sqNorm() * meas(G)));
528
529 return IFaceErr;
530 }
531}
532
533template <class T>
534void gsBiharmonicExprAssembler<T>::setSpaceBasis(const gsFunctionSet<T> & spaceBasis)
535{
536 m_spaceBasis = &spaceBasis;
537 this->_getOptions();
538 this->_initialize();
539}
540
541template <class T>
542void gsBiharmonicExprAssembler<T>::_setMapperForBiharmonic( const gsBoundaryConditions<T> & bc,
543 const gsFunctionSet<T> & spaceBasis,
544 gsDofMapper & mapper)
545{
546 if (const gsMappedBasis<2,T> * bb2 = dynamic_cast<const gsMappedBasis<2,T> *>(&spaceBasis))
547 {
548 mapper.setIdentity(bb2->nPatches(), bb2->size(), 1);
549
550 gsMatrix<index_t> bnd;
551 for (typename gsBoundaryConditions<T>::const_iterator it = bc.begin("Dirichlet"); it != bc.end("Dirichlet"); ++it)
552 {
553 bnd = bb2->basis(it->ps.patch).boundary(it->ps.side());
554 mapper.markBoundary(it->ps.patch, bnd, 0);
555 }
556
557 for (typename gsBoundaryConditions<T>::const_iterator it = bc.begin("Neumann"); it != bc.end("Neumann"); ++it)
558 {
559 bnd = bb2->basis(it->ps.patch).boundaryOffset(it->ps.side(),1);
560 mapper.markBoundary(it->ps.patch, bnd, 0);
561 }
562 mapper.finalize();
563 }
564 else if (const gsMultiBasis<T> * dbasis = dynamic_cast<const gsMultiBasis<T> *>(&spaceBasis))
565 {
566 mapper.init(*dbasis);
567
568 for (gsBoxTopology::const_iiterator it = dbasis->topology().iBegin();
569 it != dbasis->topology().iEnd(); ++it) // C^0 at the interface
570 {
571 dbasis->matchInterface(*it, mapper);
572 }
573
574 gsMatrix<index_t> bnd;
575 for (typename gsBoundaryConditions<T>::const_iterator
576 it = bc.begin("Dirichlet"); it != bc.end("Dirichlet"); ++it)
577 {
578 bnd = dbasis->basis(it->ps.patch).boundary(it->ps.side());
579 mapper.markBoundary(it->ps.patch, bnd, 0);
580 }
581
582 for (typename gsBoundaryConditions<T>::const_iterator
583 it = bc.begin("Neumann"); it != bc.end("Neumann"); ++it)
584 {
585 bnd = dbasis->basis(it->ps.patch).boundaryOffset(it->ps.side(),1);
586 mapper.markBoundary(it->ps.patch, bnd, 0);
587 }
588 mapper.finalize();
589 }
590 else
591 {
592 GISMO_ERROR("bb2 is not a gsMappedBasis");
593 }
594};
595
596
597template <class T>
598void gsBiharmonicExprAssembler<T>::_getDirichletNeumannValuesL2Projection(
599 const gsMultiPatch<T> & mp,
600 const gsMultiBasis<T> & dbasis,
601 const gsBoundaryConditions<T> & bc,
602 const gsFunctionSet<T> & spaceBasis,
603 const expr::gsFeSpace<T> & u
604 )
605{
606 if (bc.dirichletSides().size()==0 && bc.neumannSides().size()==0)
607 return;
608 if (const gsMappedBasis<2,T> * bb2 = dynamic_cast<const gsMappedBasis<2,T> *>(&spaceBasis))
609 {
610 const gsDofMapper & mapper = u.mapper();
611
612 gsMatrix<index_t> bnd = mapper.findFree(mapper.numPatches()-1);
613 gsDofMapper mapperBdy;
614 mapperBdy.setIdentity(bb2->nPatches(), bb2->size(), 1); // bb2->nPatches() == 1
615 mapperBdy.markBoundary(0, bnd, 0);
616 mapperBdy.finalize();
617
618 gsExprAssembler<T> A(1,1);
619 A.setIntegrationElements(dbasis);
620
621 auto G = A.getMap(mp);
622 auto uu = A.getSpace(*bb2);
623 auto g_bdy = A.getBdrFunction(G);
624
625 uu.setupMapper(mapperBdy);
626 gsMatrix<T> & fixedDofs_A = const_cast<expr::gsFeSpace<T>&>(uu).fixedPart();
627 fixedDofs_A.setZero( uu.mapper().boundarySize(), 1 );
628
629 A.initSystem();
630 A.assembleBdr(bc.get("Dirichlet"), uu * uu.tr() * meas(G));
631 A.assembleBdr(bc.get("Dirichlet"), uu * g_bdy * meas(G));
632 A.assembleBdr(bc.get("Neumann"),m_lambda * (igrad(uu, G) * nv(G).normalized()) * (igrad(uu, G) * nv(G).normalized()).tr() * meas(G));
633 A.assembleBdr(bc.get("Neumann"),m_lambda * (igrad(uu, G) * nv(G).normalized()) * (g_bdy.tr() * nv(G).normalized()) * meas(G));
634
635 typename gsSparseSolver<T>::SimplicialLDLT solver;
636 solver.compute( A.matrix() );
637 gsMatrix<T> & fixedDofs = const_cast<expr::gsFeSpace<T>& >(u).fixedPart();
638 fixedDofs = solver.solve(A.rhs());
639 }
640 else if (dynamic_cast<const gsMultiBasis<T> *>(&spaceBasis)) // assumes spacebasis==dbasis
641 {
642 gsDofMapper mapper = u.mapper();
643 gsDofMapper mapperBdy(dbasis, u.dim());
644 for (gsBoxTopology::const_iiterator it = dbasis.topology().iBegin();
645 it != dbasis.topology().iEnd(); ++it) // C^0 at the interface
646 {
647 dbasis.matchInterface(*it, mapperBdy);
648 }
649 for (index_t np = 0; np < mp.nPieces(); np++)
650 {
651 gsMatrix<index_t> bnd = mapper.findFree(np);
652 mapperBdy.markBoundary(np, bnd, 0);
653 }
654 mapperBdy.finalize();
655
656 gsExprAssembler<T> A(1,1);
657 A.setIntegrationElements(dbasis);
658
659 auto G = A.getMap(mp);
660 auto uu = A.getSpace(dbasis);
661 auto g_bdy = A.getBdrFunction(G);
662
663 uu.setupMapper(mapperBdy);
664 gsMatrix<T> & fixedDofs_A = const_cast<expr::gsFeSpace<T>&>(uu).fixedPart();
665 fixedDofs_A.setZero( uu.mapper().boundarySize(), 1 );
666
667 A.initSystem();
668 A.assembleBdr(bc.get("Dirichlet"), uu * uu.tr() * meas(G));
669 A.assembleBdr(bc.get("Dirichlet"), uu * g_bdy * meas(G));
670 A.assembleBdr(bc.get("Neumann"),
671 m_lambda * (igrad(uu, G) * nv(G).normalized()) * (igrad(uu, G) * nv(G).normalized()).tr() * meas(G));
672 A.assembleBdr(bc.get("Neumann"),
673 m_lambda * (igrad(uu, G) * nv(G).normalized()) * (g_bdy.tr() * nv(G).normalized()) * meas(G));
674
675 typename gsSparseSolver<T>::SimplicialLDLT solver;
676 solver.compute( A.matrix() );
677 gsMatrix<T> & fixedDofs = const_cast<expr::gsFeSpace<T>& >(u).fixedPart();
678 gsMatrix<T> fixedDofs_temp = solver.solve(A.rhs());
679
680 // Reordering the dofs of the boundary
681 fixedDofs.setZero(mapper.boundarySize(),1);
682 index_t sz = 0;
683 for (index_t np = 0; np < mp.nPieces(); np++)
684 {
685 gsMatrix<index_t> bnd = mapperBdy.findFree(np);
686 bnd.array() += sz;
687 for (index_t i = 0; i < bnd.rows(); i++)
688 {
689 index_t ii = mapperBdy.asVector()(bnd(i,0));
690 fixedDofs(mapper.global_to_bindex(mapper.asVector()(bnd(i,0))),0) = fixedDofs_temp(ii,0);
691 }
692 sz += mapperBdy.patchSize(np,0);
693 }
694 }
695 else
696 {
697 GISMO_ERROR("bb2 is not a gsMappedBasis");
698 }
699};
700
701template <class T>
702void gsBiharmonicExprAssembler<T>::_computeStabilityParameter(
703 const gsMultiPatch<T> & mp,
704 const gsMultiBasis<T> & dbasis,
705 gsMatrix<T> & mu_interfaces)
706{
707 mu_interfaces.setZero(mp.nInterfaces(),1);
708
709 index_t i = 0;
710 for ( typename gsMultiPatch<T>::const_iiterator it = mp.iBegin(); it != mp.iEnd(); ++it, ++i)
711 {
712 gsMultiPatch<T> mp_temp;
713 mp_temp.addPatch(mp.patch(it->first().patch));
714 mp_temp.addPatch(mp.patch(it->second().patch));
715 mp_temp.computeTopology();
716
717 gsMultiBasis<T> dbasis_temp;
718 dbasis_temp.addBasis(dbasis.basis(it->first().patch).clone().release());
719 dbasis_temp.addBasis(dbasis.basis(it->second().patch).clone().release());
720
721 gsBoundaryConditions<T> bc;
722
723// patchSide pS1 = mp_temp.interfaces()[0].first();
724// patchSide pS2 = mp_temp.interfaces()[0].second();
725//
726//
727// index_t side = pS1.index() < 3 ? (pS1.index() == 1 ? 2 : 1) : (pS1.index() == 3 ? 4 : 3);
728// bc.addCondition(patchSide(pS1.patchIndex(), side), condition_type::dirichlet, 0);
729//
730// side = pS2.index() < 3 ? (pS2.index() == 1 ? 2 : 1) : (pS2.index() == 3 ? 4 : 3);
731// bc.addCondition(patchSide(pS2.patchIndex(), side), condition_type::dirichlet, 0);
732
733 // Make the Eigenvalue problem to a homogeneous one
734 for (typename gsMultiPatch<T>::const_biterator bit = mp_temp.bBegin(); bit != mp_temp.bEnd(); ++bit)
735 bc.addCondition(*bit, condition_type::dirichlet, 0);
736
737 gsExprAssembler<T> A2(1, 1), B2(1, 1);
738
739 // Elements used for numerical integration
740 A2.setIntegrationElements(dbasis_temp);
741 B2.setIntegrationElements(dbasis_temp);
742
743 // Set the geometry map
744 auto GA = A2.getMap(mp_temp);
745 auto GB = B2.getMap(mp_temp);
746
747 // Set the discretization space
748 auto uA = A2.getSpace(dbasis_temp);
749 auto uB = B2.getSpace(dbasis_temp);
750
751 uA.setup(bc, dirichlet::homogeneous, 0);
752 uB.setup(bc, dirichlet::homogeneous,0);
753 //uA.setup(0);
754 //uB.setup(0);
755
756 A2.initSystem();
757 B2.initSystem();
758
759 T c = 0.25;
760 A2.assembleIfc(mp_temp.interfaces(),
761 c * ilapl(uA.left(), GA.left()) * ilapl(uA.left(), GA.left()).tr() * nv(GA.left()).norm(),
762 c * ilapl(uA.left(), GA.left()) * ilapl(uA.right(), GA.right()).tr() * nv(GA.left()).norm(),
763 c * ilapl(uA.right(), GA.right()) * ilapl(uA.left(), GA.left()).tr() * nv(GA.left()).norm(),
764 c * ilapl(uA.right(), GA.right()) * ilapl(uA.right(), GA.right()).tr() * nv(GA.left()).norm());
765
766 B2.assemble(ilapl(uB, GB) * ilapl(uB, GB).tr() * meas(GB));
767
768 // TODO INSTABLE && SLOW
769 gsMatrix<T> AA = A2.matrix().toDense();
770 gsMatrix<T> BB = B2.matrix().toDense();
771
772#ifdef gsSpectra_ENABLED
773 gsSpectraGenSymSolver<gsSparseMatrix<T>,Spectra::GEigsMode::Cholesky> ges(A2.matrix(),B2.matrix(),1,10);
774 ges.compute(Spectra::SortRule::LargestMagn,1000,1e-6,Spectra::SortRule::LargestMagn);
775 T maxEV = ges.eigenvalues()(0,0);
776#else
777 gsEigen::GeneralizedSelfAdjointEigenSolver<typename gsMatrix<T>::Base > ges(AA, BB);
778 T maxEV = ges.eigenvalues().array().maxCoeff();
779#endif
780 T m_h = dbasis_temp.basis(0).getMinCellLength(); // dbasis.basis(0).getMinCellLength();
781 mu_interfaces(i,0) = 16.0 * m_h * maxEV;
782/*
783 gsSparseSolver<T>::SimplicialLDLT sol;
784 sol.compute(B2.matrix());
785 gsSparseMatrix<T> R = sol.matrixU();
786 gsSparseMatrix<T> RT = sol.matrixL();
787 gsMatrix<T> AAA = RT.toDense().inverse() * AA * R.toDense().inverse();
788
789 gsConjugateGradient<T> cg(AAA);
790
791 cg.setCalcEigenvalues(true);
792 cg.setTolerance(1e-15);
793 cg.setMaxIterations(100000);
794
795 gsMatrix<T> rhs, result;
796 rhs.setRandom( AAA.rows(), 1 );
797 result.setRandom( AAA.rows(), 1 );
798
799 cg.solve(rhs,result);
800
801 gsInfo << "Tol: " << cg.error() << "\n";
802 gsInfo << "Max it: " << cg.iterations() << "\n";
803
804 gsMatrix<T> eigenvalues;
805 cg.getEigenvalues(eigenvalues);
806
807 gsInfo << "Cond Number: " << eigenvalues.bottomRows(1)(0,0) << "\n";
808*/
809 }
810}
811
812
813}// namespace gismo
Provides declaration of Basis abstract interface.
Provides assembler for a (planar) Biharmonic equation.
Provides gsBoundaryConditions class.
#define index_t
Definition gsConfig.h:32
#define gsDebug
Definition gsDebug.h:61
#define GISMO_ERROR(message)
Definition gsDebug.h:118
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Provides declaration of FunctionExpr class.
Provides declaration of Basis abstract interface.
Provides declaration of Basis abstract interface.
Provides declaration of a gsPiecewiseFunction class.
Header file for using Spectra extension.
EIGEN_STRONG_INLINE meas_expr< T > meas(const gsGeometryMap< T > &G)
The measure of a geometry map.
Definition gsExpressions.h:4555
EIGEN_STRONG_INLINE onormal_expr< T > nv(const gsGeometryMap< T > &u)
The (outer pointing) boundary normal of a geometry map.
Definition gsExpressions.h:4506
The G+Smo namespace, containing all definitions for the library.
S give(S &x)
Definition gsMemory.h:266