Nektar++
LinearElasticSystem.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LinearElasticSystem.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: LinearElasticSystem solve routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <algorithm>
36#include <iomanip>
37
41#include <LocalRegions/SegExp.h>
47
48#ifdef NEKTAR_USE_MPI
50#endif
51
52#ifdef NEKTAR_USE_PETSC
54#endif
55
56using namespace std;
57
58namespace Nektar
59{
60
63 "LinearElasticSystem", LinearElasticSystem::create);
64
65/*
66 * @brief Generate mapping between a simplex and the reference element.
67 *
68 * Mapping that requires the coordinates of the three vertices (x,y) of the
69 * triangular element and outputs the function coefficients that defines the
70 * mapping to the reference element for each coordinate (xfunc and yfunc)
71 *
72 * @param x co-ordinates of triangle/tetrahedron
73 * @param xf components of mapping
74 */
76{
77 int n = geom->GetNumVerts(), i, j;
78
79 DNekMat map(n, n, 1.0, eFULL);
80 DNekMat mapref(n, n, 1.0, eFULL);
81
82 // Extract coordinate information.
83 for (i = 0; i < n; ++i)
84 {
85 for (j = 0; j < n - 1; ++j)
86 {
87 map(j, i) = (*geom->GetVertex(i))[j];
88 }
89 }
90
91 // Set up reference triangle or tetrahedron mapping.
92 if (n == 3)
93 {
94 mapref(0, 0) = -1.0;
95 mapref(1, 0) = -1.0;
96 mapref(0, 1) = 1.0;
97 mapref(1, 1) = -1.0;
98 mapref(0, 2) = -1.0;
99 mapref(1, 2) = 1.0;
100 }
101 else if (n == 4)
102 {
103 mapref(0, 0) = -1.0;
104 mapref(1, 0) = -1.0;
105 mapref(2, 0) = -1.0;
106 mapref(0, 1) = 1.0;
107 mapref(1, 1) = -1.0;
108 mapref(2, 1) = -1.0;
109 mapref(0, 2) = -1.0;
110 mapref(1, 2) = 1.0;
111 mapref(2, 2) = -1.0;
112 mapref(0, 3) = -1.0;
113 mapref(1, 3) = -1.0;
114 mapref(2, 3) = 1.0;
115 }
116
117 map.Invert();
118
119 DNekMat newmap = mapref * map;
120 DNekMat mapred(n - 1, n - 1, 1.0, eFULL);
121
122 for (i = 0; i < n - 1; ++i)
123 {
124 for (j = 0; j < n - 1; ++j)
125 {
126 mapred(i, j) = newmap(i, j);
127 }
128 }
129
130 return mapred;
131}
132
133/**
134 * @brief Default constructor.
135 */
139 : EquationSystem(pSession, pGraph)
140{
141}
142
143/**
144 * @brief Set up the linear elasticity system.
145 *
146 * This routine loads the E and nu variables from the session file, creates a
147 * coupled assembly map and creates containers for the statically condensed
148 * block matrix system.
149 */
150void LinearElasticSystem::v_InitObject(bool DeclareFields)
151{
152 EquationSystem::v_InitObject(DeclareFields);
153
154 const int nVel = m_fields[0]->GetCoordim(0);
155 int n;
156
157 ASSERTL0(nVel > 1, "Linear elastic solver not set up for"
158 " this dimension (only 2D/3D supported).");
159
160 // Make sure that we have Young's modulus and Poisson ratio set.
161 m_session->LoadParameter("E", m_E, 1.00);
162 m_session->LoadParameter("nu", m_nu, 0.25);
163 m_session->LoadParameter("Beta", m_beta, 1.00);
164
165 // Create a coupled assembly map which allows us to tie u, v and w
166 // fields together.
168 std::dynamic_pointer_cast<MultiRegions::ContField>(m_fields[0]);
170 m_session, m_graph, u->GetLocalToGlobalMap(),
171 m_fields[0]->GetBndConditions(), m_fields);
172
173 // Figure out size of our new matrix systems by looping over all expansions
174 // and multiply number of coefficients by velocity components.
175 const int nEl = m_fields[0]->GetExpSize();
177
178 Array<OneD, unsigned int> sizeBnd(nEl);
179 Array<OneD, unsigned int> sizeInt(nEl);
180
181 // Allocate storage for boundary and interior map caches.
184
185 // Cache interior and boundary maps to avoid recalculating
186 // constantly. Should really be handled by manager in StdRegions but not
187 // really working yet.
188 for (n = 0; n < nEl; ++n)
189 {
190 exp = m_fields[0]->GetExp(n);
191 sizeBnd[n] = nVel * exp->NumBndryCoeffs();
192 sizeInt[n] = nVel * exp->GetNcoeffs() - sizeBnd[n];
193 exp->GetBoundaryMap(m_bmap[n]);
194 exp->GetInteriorMap(m_imap[n]);
195 }
196
197 // Create block matrix storage for the statically condensed system.
198 MatrixStorage blkmatStorage = eDIAGONAL;
200 sizeBnd, sizeBnd, blkmatStorage);
202 blkmatStorage);
204 blkmatStorage);
206 blkmatStorage);
207}
208
209/**
210 * @brief Build matrix system for linear elasticity equations.
211 *
212 * This routine constructs the matrix discretisation arising from the weak
213 * formulation of the linear elasticity equations. The resulting matrix system
214 * is then passed to LinearElasticSystem::SetStaticCondBlock in order to
215 * construct the statically condensed system.
216 *
217 * All of the matrices involved in the construction of the divergence of the
218 * stress tensor are Laplacian-like. We use the variable coefficient
219 * functionality within the library to multiply the Laplacian by the appropriate
220 * constants for all diagonal terms, and off-diagonal terms use
221 * LinearElasticSystem::BuildLaplacianIJMatrix to construct matrices containing
222 *
223 * \f[ \int \partial_{x_i} \phi_k \partial_{x_j} \phi_l \f]
224 *
225 * where mixed derivative terms are present. Symmetry (in terms of k,l) is
226 * exploited to avoid constructing this matrix repeatedly.
227 *
228 * @todo Make static condensation optional and construct full system instead.
229 */
231{
232 const int nEl = m_fields[0]->GetExpSize();
233 const int nVel = m_fields[0]->GetCoordim(0);
234
236 int n;
237
238 // Factors map for matrix keys.
240
241 // Calculate various constants
242 NekDouble mu = m_E * 0.5 / (1.0 + m_nu);
243 NekDouble lambda = m_E * m_nu / (1.0 + m_nu) / (1.0 - 2.0 * m_nu);
244
245 bool verbose = m_session->DefinesCmdLineArgument("verbose");
246 bool root = m_session->GetComm()->GetRank() == 0;
247
248 // Loop over each element and construct matrices.
249 if (nVel == 2)
250 {
251 for (n = 0; n < nEl; ++n)
252 {
253 exp = m_fields[0]->GetExp(n);
254 const int nPhys = exp->GetTotPoints();
255 Array<OneD, NekDouble> aArr(nPhys, lambda + 2 * mu);
256 Array<OneD, NekDouble> bArr(nPhys, mu);
257
258 StdRegions::VarCoeffMap varcoeffA, varcoeffD;
259 varcoeffA[StdRegions::eVarCoeffD00] = aArr;
260 varcoeffA[StdRegions::eVarCoeffD11] = bArr;
261 varcoeffD[StdRegions::eVarCoeffD00] = bArr;
262 varcoeffD[StdRegions::eVarCoeffD11] = aArr;
263
265 exp->DetShapeType(), *exp, factors,
266 varcoeffA);
268 exp->DetShapeType(), *exp, factors,
269 varcoeffD);
270
271 /*
272 * mat holds the linear operator [ A B ] acting on [ u ].
273 * [ C D ] [ v ]
274 */
276 mat[0][0] = exp->GenMatrix(matkeyA);
277 mat[0][1] = BuildLaplacianIJMatrix(1, 0, mu + lambda, exp);
278 mat[1][0] = mat[0][1];
279 mat[1][1] = exp->GenMatrix(matkeyD);
280
281 // Set up the statically condensed block for this element.
282 SetStaticCondBlock(n, exp, mat);
283
284 if (verbose && root)
285 {
286 cout << "\rBuilding matrix system: " << (int)(100.0 * n / nEl)
287 << "%" << flush;
288 }
289 }
290 }
291 else if (nVel == 3)
292 {
293 for (n = 0; n < nEl; ++n)
294 {
295 exp = m_fields[0]->GetExp(n);
296 const int nPhys = exp->GetTotPoints();
297 Array<OneD, NekDouble> aArr(nPhys, lambda + 2 * mu);
298 Array<OneD, NekDouble> bArr(nPhys, mu);
299
300 StdRegions::VarCoeffMap varcoeffA, varcoeffE, varcoeffI;
301 varcoeffA[StdRegions::eVarCoeffD00] = aArr;
302 varcoeffA[StdRegions::eVarCoeffD11] = bArr;
303 varcoeffA[StdRegions::eVarCoeffD22] = bArr;
304 varcoeffE[StdRegions::eVarCoeffD00] = bArr;
305 varcoeffE[StdRegions::eVarCoeffD11] = aArr;
306 varcoeffE[StdRegions::eVarCoeffD22] = bArr;
307 varcoeffI[StdRegions::eVarCoeffD00] = bArr;
308 varcoeffI[StdRegions::eVarCoeffD11] = bArr;
309 varcoeffI[StdRegions::eVarCoeffD22] = aArr;
310
312 exp->DetShapeType(), *exp, factors,
313 varcoeffA);
315 exp->DetShapeType(), *exp, factors,
316 varcoeffE);
318 exp->DetShapeType(), *exp, factors,
319 varcoeffI);
320
321 /*
322 * mat holds the linear operator [ A B C ] acting on [ u ].
323 * [ D E F ] [ v ]
324 * [ G H I ] [ w ]
325 */
327 mat[0][0] = exp->GenMatrix(matkeyA);
328 mat[0][1] = BuildLaplacianIJMatrix(1, 0, mu + lambda, exp);
329 mat[0][2] = BuildLaplacianIJMatrix(2, 0, mu + lambda, exp);
330
331 mat[1][0] = mat[0][1];
332 mat[1][1] = exp->GenMatrix(matkeyE);
333 mat[1][2] = BuildLaplacianIJMatrix(2, 1, mu + lambda, exp);
334
335 mat[2][0] = mat[0][2];
336 mat[2][1] = mat[1][2];
337 mat[2][2] = exp->GenMatrix(matkeyI);
338
339 // Set up the statically condensed block for this element.
340 SetStaticCondBlock(n, exp, mat);
341
342 if (verbose && root)
343 {
344 cout << "\rBuilding matrix system: " << (int)(100.0 * n / nEl)
345 << "%" << flush;
346 }
347 }
348 }
349
350 if (verbose && root)
351 {
352 cout << "\rBuilding matrix system: done." << endl;
353 }
354}
355
356/**
357 * @brief Generate summary at runtime.
358 */
360{
362
363 AddSummaryItem(s, "Young's modulus", m_E);
364 AddSummaryItem(s, "Poisson ratio", m_nu);
365}
366
367/**
368 * @brief Solve elliptic linear elastic system.
369 *
370 * The solve proceeds as follows:
371 *
372 * - Create a MultiRegions::GlobalLinSys object.
373 * - Evaluate a forcing function from the session file.
374 * - If the temperature term is enabled, evaluate this and add to forcing.
375 * - Apply Dirichlet boundary conditions.
376 * - Scatter forcing into the correct ordering according to the coupled
377 * assembly map.
378 * - Do the solve.
379 * - Scatter solution back to fields and backwards transform to physical space.
380 */
382{
383 int i, j, k, l, nv;
384 const int nVel = m_fields[0]->GetCoordim(0);
385
386 // Build initial matrix system.
388
389 // Now we've got the matrix system set up, create a GlobalLinSys
390 // object. We mask ourselves as LinearAdvectionReaction to create a full
391 // matrix instead of symmetric storage.
395
396 // Currently either direct or iterative static condensation is
397 // supported.
398 if (m_assemblyMap->GetGlobalSysSolnType() ==
400 {
404 }
405 else if (m_assemblyMap->GetGlobalSysSolnType() ==
407 {
412 }
413#ifdef NEKTAR_USE_PETSC
414 else if (m_assemblyMap->GetGlobalSysSolnType() ==
416 {
420 }
421#endif
422#ifdef NEKTAR_USE_MPI
423 else if (m_assemblyMap->GetGlobalSysSolnType() ==
425 {
429 }
430#endif
431
432 linSys->Initialise(m_assemblyMap);
433
434 const int nCoeffs = m_fields[0]->GetNcoeffs();
435
436 //
437 // -- Evaluate forcing functions
438 //
439
440 // Evaluate the forcing function from the XML file.
442 GetFunction("Forcing")->Evaluate(forcing);
443
444 // Add temperature term
445 string tempEval;
446 m_session->LoadSolverInfo("Temperature", tempEval, "None");
447
448 if (tempEval == "Jacobian")
449 {
450 // Allocate storage
452
453 for (nv = 0; nv < nVel; ++nv)
454 {
456 m_temperature[nv] =
458
459 for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
460 {
461 // Calculate element area
462 LocalRegions::ExpansionSharedPtr exp = m_fields[0]->GetExp(i);
463 LibUtilities::PointsKeyVector pkey = exp->GetPointsKeys();
464 Array<OneD, NekDouble> jac = exp->GetMetricInfo()->GetJac(pkey);
465
466 int offset = m_fields[0]->GetPhys_Offset(i);
467
468 if (exp->GetMetricInfo()->GetGtype() ==
470 {
471 Vmath::Smul(exp->GetTotPoints(), m_beta, jac, 1,
472 tmp = m_temperature[nv] + offset, 1);
473 }
474 else
475 {
476 Vmath::Fill(exp->GetTotPoints(), m_beta * jac[0],
477 tmp = m_temperature[nv] + offset, 1);
478 }
479 }
480 m_fields[nv]->PhysDeriv(nv, m_temperature[nv], forcing[nv]);
481 }
482 }
483 else if (tempEval == "Metric")
484 {
485 ASSERTL0((m_fields[0]->GetCoordim(0) == 2 &&
486 m_graph->GetAllQuadGeoms().size() == 0) ||
487 (m_fields[0]->GetCoordim(0) == 3 &&
488 m_graph->GetAllPrismGeoms().size() == 0 &&
489 m_graph->GetAllPyrGeoms().size() == 0 &&
490 m_graph->GetAllHexGeoms().size() == 0),
491 "LinearIdealMetric temperature only implemented for "
492 "two-dimensional triangular meshes or three-dimensional "
493 "tetrahedral meshes.");
494
497
498 for (nv = 0; nv < nVel; ++nv)
499 {
500 m_temperature[nv] =
503
504 for (i = 0; i < nVel; ++i)
505 {
506 m_stress[nv][i] =
508 }
509 }
510
511 for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
512 {
513 // Calculate element area
514 LocalRegions::ExpansionSharedPtr exp = m_fields[0]->GetExp(i);
515 LibUtilities::PointsKeyVector pkey = exp->GetPointsKeys();
517 exp->GetMetricInfo()->GetDeriv(pkey);
518 int offset = m_fields[0]->GetPhys_Offset(i);
519
520 DNekMat i2rm = MappingIdealToRef(exp->GetGeom());
521
522 // Compute metric tensor
523 DNekMat jac(nVel, nVel, 0.0, eFULL);
524 DNekMat jacIdeal(nVel, nVel, 0.0, eFULL);
525 DNekMat metric(nVel, nVel, 0.0, eFULL);
526
527 for (j = 0; j < deriv[0][0].size(); ++j)
528 {
529 for (k = 0; k < nVel; ++k)
530 {
531 for (l = 0; l < nVel; ++l)
532 {
533 jac(l, k) = deriv[k][l][j];
534 }
535 }
536
537 jacIdeal = jac * i2rm;
538 metric = Transpose(jacIdeal) * jacIdeal;
539
540 // Compute eigenvalues/eigenvectors of metric tensor using
541 // ideal mapping.
542 char jobvl = 'N', jobvr = 'V';
543 int worklen = 8 * nVel, info;
544
545 DNekMat eval(nVel, nVel, 0.0, eDIAGONAL);
546 DNekMat evec(nVel, nVel, 0.0, eFULL);
547 DNekMat evecinv(nVel, nVel, 0.0, eFULL);
548 Array<OneD, NekDouble> vl(nVel * nVel);
549 Array<OneD, NekDouble> work(worklen);
550 Array<OneD, NekDouble> wi(nVel);
551
552 Lapack::Dgeev(jobvl, jobvr, nVel, metric.GetRawPtr(), nVel,
553 &(eval.GetPtr())[0], &wi[0], &vl[0], nVel,
554 &(evec.GetPtr())[0], nVel, &work[0], worklen,
555 info);
556
557 evecinv = evec;
558 evecinv.Invert();
559
560 // rescaling of the eigenvalues
561 for (nv = 0; nv < nVel; ++nv)
562 {
563 eval(nv, nv) = m_beta * (sqrt(eval(nv, nv)) - 1.0);
564 }
565
566 DNekMat beta = evec * eval * evecinv;
567
568 for (k = 0; k < nVel; ++k)
569 {
570 for (l = 0; l < nVel; ++l)
571 {
572 m_stress[k][l][offset + j] = beta(k, l);
573 }
574 }
575 }
576
577 if (deriv[0][0].size() != exp->GetTotPoints())
578 {
580 for (k = 0; k < nVel; ++k)
581 {
582 for (l = 0; l < nVel; ++l)
583 {
584 Vmath::Fill(exp->GetTotPoints(), m_stress[k][l][offset],
585 tmp = m_stress[k][l] + offset, 1);
586 }
587 }
588 }
589 }
590
591 // Calculate divergence of stress tensor.
593 for (i = 0; i < nVel; ++i)
594 {
595 for (j = 0; j < nVel; ++j)
596 {
597 m_fields[i]->PhysDeriv(j, m_stress[i][j], tmpderiv);
598 Vmath::Vadd(m_fields[i]->GetNpoints(), tmpderiv, 1,
599 m_temperature[i], 1, m_temperature[i], 1);
600 }
601
603 forcing[i], 1);
604 }
605 }
606 else if (tempEval != "None")
607 {
608 ASSERTL0(false, "Unknown temperature form: " + tempEval);
609 }
610
611 // Set up some temporary storage.
612 //
613 // - forCoeffs holds the forcing coefficients in a local ordering;
614 // however note that the ordering is different and dictated by the
615 // assembly map. We loop over each element, then the boundary degrees of
616 // freedom for u, boundary for v, followed by the interior for u and then
617 // interior for v.
618 // - inout holds the Dirichlet degrees of freedom in the local ordering,
619 // which have the boundary conditions imposed
620 Array<OneD, NekDouble> forCoeffs(nVel * nCoeffs, 0.0);
621 Array<OneD, NekDouble> inout(nVel * nCoeffs, 0.0);
622 Array<OneD, NekDouble> tmp(nCoeffs);
623
624 for (nv = 0; nv < nVel; ++nv)
625 {
626 // Take the inner product of the forcing function.
627 m_fields[nv]->IProductWRTBase(forcing[nv], tmp);
628
629 // Impose Dirichlet condition on field which should be initialised
630 Array<OneD, NekDouble> loc_inout = m_fields[nv]->UpdateCoeffs();
631 m_fields[nv]->ImposeDirichletConditions(loc_inout);
632
633 // Scatter forcing into RHS vector according to the ordering dictated in
634 // the comment above.
635 for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
636 {
637 int nBnd = m_bmap[i].size();
638 int nInt = m_imap[i].size();
639 int offset = m_fields[nv]->GetCoeff_Offset(i);
640
641 for (j = 0; j < nBnd; ++j)
642 {
643 forCoeffs[nVel * offset + nv * nBnd + j] =
644 -1.0 * tmp[offset + m_bmap[i][j]];
645 inout[nVel * offset + nv * nBnd + j] =
646 loc_inout[offset + m_bmap[i][j]];
647 }
648 for (j = 0; j < nInt; ++j)
649 {
650 forCoeffs[nVel * (offset + nBnd) + nv * nInt + j] =
651 -1.0 * tmp[offset + m_imap[i][j]];
652 }
653 }
654 }
655
656 //
657 // -- Perform solve
658 //
659 // Solve.
660 linSys->Solve(forCoeffs, inout, m_assemblyMap);
661
662 //
663 // -- Postprocess
664 //
665
666 // Scatter back to field degrees of freedom
667 for (nv = 0; nv < nVel; ++nv)
668 {
669 for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
670 {
671 int nBnd = m_bmap[i].size();
672 int nInt = m_imap[i].size();
673 int offset = m_fields[nv]->GetCoeff_Offset(i);
674
675 for (j = 0; j < nBnd; ++j)
676 {
677 m_fields[nv]->UpdateCoeffs()[offset + m_bmap[i][j]] =
678 inout[nVel * offset + nv * nBnd + j];
679 }
680 for (j = 0; j < nInt; ++j)
681 {
682 m_fields[nv]->UpdateCoeffs()[offset + m_imap[i][j]] =
683 inout[nVel * (offset + nBnd) + nv * nInt + j];
684 }
685 }
686 m_fields[nv]->BwdTrans(m_fields[nv]->GetCoeffs(),
687 m_fields[nv]->UpdatePhys());
688 }
689}
690
691/**
692 * @brief Given a block matrix for an element, construct its static condensation
693 * matrices.
694 *
695 * This routine essentially duplicates the logic present in many of the
696 * LocalRegions matrix generation routines to construct the statically condensed
697 * equivalent of mat to pass to the GlobalLinSys solver.
698 *
699 * @param n Element/block number.
700 * @param exp Pointer to expansion.
701 * @param mat Block matrix containing matrix operator.
702 */
704 const int n, const LocalRegions::ExpansionSharedPtr exp,
706{
707 int i, j, k, l;
708 const int nVel = mat.GetRows();
709 const int nB = exp->NumBndryCoeffs();
710 const int nI = exp->GetNcoeffs() - nB;
711 const int nBnd = exp->NumBndryCoeffs() * nVel;
712 const int nInt = exp->GetNcoeffs() * nVel - nBnd;
713 const MatrixStorage s = eFULL; // Maybe look into doing symmetric
714 // version of this?
715
724
725 for (i = 0; i < nVel; ++i)
726 {
727 for (j = 0; j < nVel; ++j)
728 {
729 // Boundary-boundary and boundary-interior
730 for (k = 0; k < nB; ++k)
731 {
732 for (l = 0; l < nB; ++l)
733 {
734 (*A)(k + i * nB, l + j * nB) =
735 (*mat[i][j])(m_bmap[n][k], m_bmap[n][l]);
736 }
737
738 for (l = 0; l < nI; ++l)
739 {
740 (*B)(k + i * nB, l + j * nI) =
741 (*mat[i][j])(m_bmap[n][k], m_imap[n][l]);
742 }
743 }
744
745 // Interior-boundary / interior-interior
746 for (k = 0; k < nI; ++k)
747 {
748 for (l = 0; l < nB; ++l)
749 {
750 (*C)(k + i * nI, l + j * nB) =
751 (*mat[i][j])(m_imap[n][k], m_bmap[n][l]);
752 }
753
754 for (l = 0; l < nI; ++l)
755 {
756 (*D)(k + i * nI, l + j * nI) =
757 (*mat[i][j])(m_imap[n][k], m_imap[n][l]);
758 }
759 }
760 }
761 }
762
763 // Construct static condensation matrices.
764 D->Invert();
765 (*B) = (*B) * (*D);
766 (*A) = (*A) - (*B) * (*C);
767
768 DNekScalMatSharedPtr tmp_mat;
769 m_schurCompl->SetBlock(
771 m_BinvD->SetBlock(
772 n, n, tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, B));
773 m_C->SetBlock(
774 n, n, tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, C));
775 m_Dinv->SetBlock(
776 n, n, tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, D));
777}
778
779/**
780 * @brief Construct a LaplacianIJ matrix for a given expansion.
781 *
782 * This routine constructs matrices whose entries contain the evaluation of
783 *
784 * \f[ \partial_{\tt k1} \phi_i \partial_{\tt k2} \phi_j \,dx \f]
785 */
787 const int k1, const int k2, const NekDouble scale,
789{
790 const int nCoeffs = exp->GetNcoeffs();
791 const int nPhys = exp->GetTotPoints();
792 int i;
793
794 DNekMatSharedPtr ret =
795 MemoryManager<DNekMat>::AllocateSharedPtr(nCoeffs, nCoeffs, 0.0, eFULL);
796
797 Array<OneD, NekDouble> tmp2(nPhys);
798 Array<OneD, NekDouble> tmp3(nPhys);
799
800 for (i = 0; i < nCoeffs; ++i)
801 {
802 Array<OneD, NekDouble> tmp1(nCoeffs, 0.0);
803 tmp1[i] = 1.0;
804
805 exp->BwdTrans(tmp1, tmp2);
806 exp->PhysDeriv(k1, tmp2, tmp3);
807 exp->IProductWRTDerivBase(k2, tmp3, tmp1);
808
809 Vmath::Smul(nCoeffs, scale, &tmp1[0], 1,
810 &(ret->GetPtr())[0] + i * nCoeffs, 1);
811 }
812
813 return ret;
814}
815
817 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
818 std::vector<std::string> &variables)
819{
820 const int nVel = m_fields[0]->GetCoordim(0);
821 const int nCoeffs = m_fields[0]->GetNcoeffs();
822
823 if (m_temperature.size() == 0)
824 {
825 return;
826 }
827
828 for (int i = 0; i < nVel; ++i)
829 {
830 Array<OneD, NekDouble> tFwd(nCoeffs);
831 m_fields[i]->FwdTrans(m_temperature[i], tFwd);
832 fieldcoeffs.push_back(tFwd);
833 variables.push_back("ThermStressDiv" + std::to_string(i));
834 }
835
836 if (m_stress.size() == 0)
837 {
838 return;
839 }
840
841 for (int i = 0; i < nVel; ++i)
842 {
843 for (int j = 0; j < nVel; ++j)
844 {
845 Array<OneD, NekDouble> tFwd(nCoeffs);
846 m_fields[i]->FwdTrans(m_stress[i][j], tFwd);
847 fieldcoeffs.push_back(tFwd);
848 variables.push_back("ThermStress" + std::to_string(i) +
849 std::to_string(j));
850 }
851 }
852}
853
854} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Array< OneD, Array< OneD, unsigned int > > m_imap
Interior maps for each of the fields.
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
DNekScalBlkMatSharedPtr m_Dinv
Matrix of elemental components.
NekDouble m_nu
Poisson ratio.
static std::string className
Name of class.
void v_InitObject(bool DeclareFields=true) override
Set up the linear elasticity system.
void BuildMatrixSystem()
Build matrix system for linear elasticity equations.
LinearElasticSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Default constructor.
void SetStaticCondBlock(const int n, const LocalRegions::ExpansionSharedPtr exp, Array< TwoD, DNekMatSharedPtr > &mat)
Given a block matrix for an element, construct its static condensation matrices.
Array< OneD, Array< OneD, NekDouble > > m_temperature
Storage for the temperature terms.
CoupledAssemblyMapSharedPtr m_assemblyMap
Assembly map for the coupled (u,v,w) system.
void v_DoSolve() override
Solve elliptic linear elastic system.
DNekScalBlkMatSharedPtr m_BinvD
Matrix of elemental components.
NekDouble m_beta
Parameter dictating amount of thermal stress to add.
DNekScalBlkMatSharedPtr m_schurCompl
Schur complement boundary-boundary matrix.
void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) override
NekDouble m_E
Young's modulus.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_stress
Storage for the thermal stress terms.
Array< OneD, Array< OneD, unsigned int > > m_bmap
Boundary maps for each of the fields.
DNekMatSharedPtr BuildLaplacianIJMatrix(const int k1, const int k2, const NekDouble scale, LocalRegions::ExpansionSharedPtr exp)
Construct a LaplacianIJ matrix for a given expansion.
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Generate summary at runtime.
DNekScalBlkMatSharedPtr m_C
Matrix of elemental components.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
A base class for describing how to solve specific equations.
SOLVER_UTILS_EXPORT void SessionSummary(SummaryList &vSummary)
Write out a session summary.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
static void Dgeev(const char &uplo, const char &lrev, const int &n, const double *a, const int &lda, double *wr, double *wi, double *rev, const int &ldr, double *lev, const int &ldv, double *work, const int &lwork, int &info)
Solve general real matrix eigenproblem.
Definition: Lapack.hpp:329
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
static PreconditionerSharedPtr NullPreconditionerSharedPtr
std::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:51
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:268
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:51
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
DNekMat MappingIdealToRef(SpatialDomains::GeometrySharedPtr geom)
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294