G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
kirchhoff-Love_example.cpp

In this example, the Kirchhoff-Love shell is implemented. References for this implementation include the original work by Kiendl et. al. (2009) [1], the PhD thesis of Josef Kiendl (2011) [2], the PhD thesis of Anmol Goyal (2015) [3], the MSc thesis of Hugo Verhelst (2019) [4]. Anmol contributed to the first Kirchhoff-Love implementation of G+Smo and Hugo to the current one.

Governing Equations

The variational formulation for the Kirchhoff-Love shell is:

\begin{eqnarray*} \int_\Omega \mathbf{f}\cdot\mathbf{v}\ + p\mathbf{\hat{n}}:\text{d}\Omega - \int_\Omega \mathbf{n}(\mathbf{u}):\mathbf{\varepsilon}^\prime(\mathbf{u},\mathbf{v}) + \mathbf{m}(\mathbf{u}):\mathbf{\kappa}^\prime(\mathbf{u},\mathbf{v}) \text{d}\Omega = 0 \end{eqnarray*}

The left-hand side of this formulation in fact is the residual of the problem, being the balance between external and internal forces. In this equation, \(\Omega\) is the domain, \( \mathbf{f} \) is a vector with a distributed load acting on the mid-plane of the shell, \( \mathbf{v} \) is the three-dimensional test function, \( p\) is a follower pressure acting on the shell, \( \mathbf{\hat{n}}\) is the shell normal, \( n \) is the normal-force tensor, \( \mathbf{\varepsilon}^\prime \) is the variation of the membrane strain tensor, \( \mathbf{m} \) is the bending moment tensor, \( \mathbf{\kappa}^\prime \) is the variation of the bending strain tensor and \(\mathbf{u} \) is the displacement of the mid-plane of the shell; thus the solution to the problem.

In order to solve the variational formulation for \( \mathbf{u} \), the second variation is derived, such that a lineararized form can be found or such that the Jacobian and for Newton iterations can be found. The second variation is:

\begin{eqnarray*} j(\mathbf{u},\mathbf{v},\mathbf{w}) = \int_\Omega \mathbf{n}(\mathbf{u}):\mathbf{\varepsilon}^{\prime\prime}(\mathbf{u},\mathbf{v},\mathbf{v}) + \mathbf{n}^\prime(\mathbf{u},\mathbf{v}):\mathbf{\varepsilon}^{\prime}(\mathbf{u},\mathbf{v}) + \mathbf{m}(\mathbf{u}):\mathbf{\kappa}^{\prime\prime}(\mathbf{u},\mathbf{v},\mathbf{w}) + \mathbf{m}^\prime(\mathbf{u},\mathbf{v}):\mathbf{\kappa}^\prime(\mathbf{u},\mathbf{v}) \text{d}\Omega - \int_\Omega p \mathbf{v}\cdot\mathbf{\hat{n}}^\prime :\text{d}\Omega \end{eqnarray*}

In this equation, the primed ( \(\cdot^\prime\)) expressions again denote the first variation and the double-primed ( \(\cdot^{\prime\prime}\)) expressions denote second variations.

From the above equations, a system of equations can be assembled, which will later be done by the gsExprAssembler . In case of a lineararized system, only the single-primed expressions are non-zero, together with the term containing the external load.

The expressions for the strains, normal foce and bending moment are as follows (see sec. 3.2 and 3.3 of [3]):

\begin{align*} \varepsilon_{\alpha\beta} &= \frac{\partial\mathbf{c}}{\partial\theta_\alpha}\cdot\frac{\partial\mathbf{c}}{\partial\theta_\beta} - \frac{\partial\mathbf{C}}{\partial\theta_\alpha}\cdot\frac{\partial\mathbf{C}}{\partial\theta_\beta}\\ \varepsilon^\prime_{\alpha\beta} &= \frac{1}{2}\frac{\partial\mathbf{v}}{\partial\theta_\alpha}\cdot\frac{\partial\mathbf{c}}{\partial\theta_\beta} + \frac{1}{2}\frac{\partial\mathbf{c}}{\partial\theta_\alpha}\cdot\frac{\partial\mathbf{v}}{\partial\theta_\beta} \\ \kappa_{\alpha\beta} &= \frac{\partial^2\mathbf{C}}{\partial\theta_\alpha\partial\theta_\beta}\cdot\mathbf{\hat{N}} - \frac{\partial^2\mathbf{c}}{\partial\theta_\alpha\partial\theta_\beta}\cdot\mathbf{\hat{n}}\\ \kappa^\prime_{\alpha\beta} &= \frac{\partial^2\mathbf{v}}{\partial\theta_\alpha\partial\theta_\beta}\cdot\mathbf{\hat{n}} + \frac{\partial^2\mathbf{c}}{\partial\theta_\alpha\partial\theta_\beta}\cdot\mathbf{\hat{n}}^\prime(\mathbf{u})\\ \mathbf{\hat{n}}^\prime(\mathbf{v}) &= \mathbf{m}_{\mathbf{v}} - (\mathbf{\hat{n}}\cdot \mathbf{m}_{\mathbf{v}})\mathbf{\hat{n}} \\ \mathbf{m}_{\mathbf{v}} &= \frac{1}{\vert \frac{\partial\mathbf{c}}{\partial\theta_1}\times \frac{\partial\mathbf{c}}{\partial\theta_2}\vert } \left( \frac{\partial\mathbf{v}}{\partial\theta_1}\times \frac{\partial\mathbf{c}}{\partial\theta_2} + \frac{\partial\mathbf{c}}{\partial\theta_1} \times \frac{\partial\mathbf{v}}{\partial\theta_2} \right)\\ \varepsilon^{\prime\prime}_{\alpha\beta} &= \frac{\partial\mathbf{v}}{\partial\theta_\alpha}\cdot\frac{\partial\mathbf{w}}{\partial\theta_\beta}\\ \kappa^{\prime\prime}_{\alpha\beta} &= \frac{\partial^2\mathbf{v}}{\partial\theta_\alpha\partial\theta_\beta}\cdot\mathbf{\hat{n}}^\prime(\mathbf{w}) + \frac{\partial^2\mathbf{w}}{\partial\theta_\alpha\partial\theta_\beta}\cdot\mathbf{\hat{n}}^\prime(\mathbf{v}) + \frac{\partial^2\mathbf{c}}{\partial\theta_\alpha\partial\theta_\beta}\cdot\mathbf{\hat{n}}^{\prime\prime}(\mathbf{v},\mathbf{w})\\ \mathbf{\hat{n}}^{\prime\prime}(\mathbf{v},\mathbf{w}) &= \mathbf{m}_{\mathbf{v}}^\prime - ( \mathbf{m}_{\mathbf{v}}\cdot\mathbf{\hat{n}}^\prime(\mathbf{w}) + \mathbf{\hat{n}}\cdot\mathbf{m}_{\mathbf{v}}^\prime(\mathbf{w}) )\mathbf{\hat{n}} - ( \mathbf{\hat{n}}\cdot \mathbf{m}_{\mathbf{v}} )\mathbf{\hat{n}}^\prime(\mathbf{w})\\ \mathbf{m}_{\mathbf{v}}^{\prime}(\mathbf{w}) &= \mathbf{m}_{\mathbf{vw}}-(\mathbf{n}\cdot\mathbf{m}_{\mathbf{w}})\mathbf{m}_{\mathbf{w}}\\ \mathbf{m}_{\mathbf{vw}} &= \frac{1}{\vert \frac{\partial\mathbf{c}}{\partial\theta_1}\times \frac{\partial\mathbf{c}}{\partial\theta_2}\vert } \left( \frac{\partial\mathbf{v}}{\partial\theta_1}\times \frac{\partial\mathbf{w}}{\partial\theta_2} + \frac{\partial\mathbf{w}}{\partial\theta_1} \times \frac{\partial\mathbf{v}}{\partial\theta_2} \right)\\ n_{\alpha\beta} &= t \mathcal{C}^{\alpha\beta\gamma\delta} \varepsilon_{\alpha\beta}\\ n_{\alpha\beta}^\prime &= \mathcal{C}^{\alpha\beta\gamma\delta} \varepsilon_{\alpha\beta}^{\prime}\\ m_{\alpha\beta} &= t \mathcal{C}^{\alpha\beta\gamma\delta} \kappa{\alpha\beta}\\ m_{\alpha\beta}^\prime &= \mathcal{C}^{\alpha\beta\gamma\delta} \kappa_{\alpha\beta}^{\prime} \end{align*}

Here, \( t\) is the thickness of the shell and \( \mathcal{C}^{\alpha\beta\gamma\delta}\) is the material tensor (see [3] eq. 3.15).

Implementation of the Kirchhoff-Love Shell

Now we will have a look at the implementation in kirchhoff-Love_example.cpp.

The file starts with a large amount of classes being defined. These classes all inherit from the class _expr which means that they are G+Smo expressions. Without going into detail for the implementations of all the classes, the following functions match with the expressions defined above. Note that all these expressions are specific for the shell problem.

  • var1_expr resembles \( \mathbf{\hat{n}}^\prime\) and requires the space and the geometry.
  • var2_expr resembles \( \mathbf{\hat{n}}^{\prime\prime} \cdot \mathbf{p} \) for an arbitrary vector \( \mathbf{p} \).
  • deriv2dot_expr resembles either \( \frac{\partial^2\mathbf{v}}{\partial\theta_\alpha\partial\theta_\beta}\cdot \mathbf{p} \) or \( \frac{\partial^2\mathbf{c}}{\partial\theta_\alpha\partial\theta_\beta}\cdot \mathbf{p} \) for arbitrary vector \( \mathbf{p} \) space \( \mathbf{v}\) and geometry \( c \),
  • deriv2_expr resembles either \( \frac{\partial^2\mathbf{v}}{\partial\theta_\alpha\partial\theta_\beta} \) or \( \frac{\partial^2\mathbf{c}}{\partial\theta_\alpha\partial\theta_\beta} \) for space \( \mathbf{v}\) or geometry \( c \),
  • flatdot_expr is used to compute \( \mathbf{\varepsilon}^{\prime\prime}\cdot\mathbf{n}\) and its main purpose is memory-efficient implementation
  • flatdot2_expr is used to compute the first two terms of \( \mathbf{\kappa}^{\prime\prime}\cdot\mathbf{m}\) and its main purpose is memory-efficient implementation
  • cartcov_expr returns a transformation matrix between local covariant to local cartesian coordinates. The expression for its inverse, cartcovinv_expr is also implemented,
  • cartcov_expr returns a transformation matrix between local contravariant to local cartesian coordinates. The expression for its inverse, cartconinv_expr is also implemented,

After the shell-specific expressions, the material matrix class is defined, which is used to handle the material tensor ${C}$. This class inherits from gsFunction because it will later be used as a function to declare a variable in the gsExprAssembler.

Now, let us skip the implementation of the shell-specific expression implementations and the definition of the material matrix class. Before the main file starts, there are some forward declarations for functions that will help to define rectangular domains with arbitrary number of elements, order and length and width.

template <class T>
gsMultiPatch<T> RectangularDomain(int n, int m, int p, int q, T L, T B);
template <class T>
gsMultiPatch<T> RectangularDomain(int n, int p, T L, T B);

The first lines of our main file will handle the command line arguments for our script. This works similar to what can be seen in poisson2_example.cpp. For this particular case, there are multiple test cases that can be chosen, the flag for nonlinear simulation can be set and the usage of weak Dirichlet boundary conditions can be triggered for some test cases.

bool plot = false;
index_t numRefine = 1;
index_t numElevate = 1;
index_t testCase = 1;
bool weak = false;
bool nonlinear = false;
std::string fn;
real_t E_modulus = 1.0;
real_t PoissonRatio = 0.0;
real_t thickness = 1.0;
gsCmdLine cmd("Tutorial on solving a Poisson problem.");
cmd.addInt( "e", "degreeElevation",
"Number of degree elevation steps to perform before solving (0: equalize degree in all directions)", numElevate );
cmd.addInt( "r", "uniformRefine", "Number of Uniform h-refinement steps to perform before solving", numRefine );
cmd.addInt( "t", "testCase", "Test case to run: 1 = unit square; 2 = Scordelis Lo Roof", testCase );
cmd.addSwitch("nl", "Solve nonlinear problem", nonlinear);
cmd.addSwitch("plot", "Create a ParaView visualization file with the solution", plot);
cmd.addSwitch("weak", "Weak BCs", weak);
try { cmd.getValues(argc,argv); } catch (int rv) { return rv; }

Next up, we define the geometry and material data for different test cases. We define two multipatches; one to store the original geometry and one to store the deformed geometry. Furthermore, each test case requires a thickness, a Youngs modulus, a Poisson's ratio and a geometry.

gsMultiPatch<> mp;
gsMultiPatch<> mp_def;
if (testCase == 2 || testCase == 3)
{
thickness = 0.25;
E_modulus = 4.32E8;
fn = "surfaces/scordelis_lo_roof.xml";
gsReadFile<>(fn, mp);
PoissonRatio = 0.0;
}
else if (testCase == 6)
{
thickness = 0.1;
real_t mu = 4.225e5;
PoissonRatio = 0.3;
E_modulus = (2+PoissonRatio)*mu;
fn = "surfaces/quarter_sphere.xml";
gsReadFile<>(fn, mp);
}
else
{
// Unit square
mp.addPatch( gsNurbsCreator<>::BSplineSquare(1) ); // degree
mp.addAutoBoundaries();
mp.embed(3);
E_modulus = 1e0;
thickness = 1e0;
PoissonRatio = 0.0;
}

In the following snippet, the geometry is refined and its degree is elevated.

// p-refine
if (numElevate!=0)
mp.degreeElevate(numElevate);
// h-refine
for (int r =0; r < numRefine; ++r)
mp.uniformRefine();
mp_def = mp;
gsWriteParaview(mp,"mp",1000,true);
gsMultiBasis<> dbasis(mp);
gsInfo << "Patches: "<< mp.nPatches() <<", degree: "<< dbasis.minCwiseDegree() <<"\n";
gsInfo << dbasis.basis(0)<<"\n";

Now we define the boundary conditions, add the geometry map to the gsBoundaryConditions class and we define vectors for the implementation of the surface load (tmp), the neumann boundary conditions (neu) and functions are added for the imposed displacements, weak boundary conditions and the Neumann boundary conditions. We also have a scalar for applied pressures.

gsBoundaryConditions<> bc;
bc.setGeoMap(mp);
gsVector<> tmp(3);
tmp << 0, 0, 0;
gsVector<> neu(3);
neu << 0, 0, 0;
gsConstantFunction<> displx(0.09317696,3);
gsConstantFunction<> neuData(neu,3);
gsConstantFunction<> weakBC(neu,3);
real_t pressure = 0.0;

Now let us consider the definition of the individual test cases. The first case is a square plate (see geometry definition) where all boundaries are fixed in all directions. Note that the boundary conditions are applied per component (last entry) of the displacement.

if (testCase == 1)
{
for (index_t i=0; i!=3; ++i)
{
bc.addCondition(boundary::north, condition_type::dirichlet, 0, 0, false, i ); // unknown 0 - x
bc.addCondition(boundary::east, condition_type::dirichlet, 0, 0, false, i ); // unknown 1 - y
bc.addCondition(boundary::south, condition_type::dirichlet, 0, 0, false, i ); // unknown 2 - z
bc.addCondition(boundary::west, condition_type::dirichlet, 0, 0, false, i ); // unknown 2 - z
}
tmp << 0,0,-1;
}

For the second test case, a corner is fixed.

else if (testCase == 2)
{
// Diaphragm conditions
bc.addCondition(boundary::west, condition_type::dirichlet, 0, 0, false, 1 ); // unknown 1 - y
bc.addCondition(boundary::west, condition_type::dirichlet, 0, 0, false, 2 ); // unknown 2 - z
bc.addCornerValue(boundary::southwest, 0.0, 0, 0, 0); // (corner,value, patch, unknown, component)
bc.addCondition(boundary::east, condition_type::dirichlet, 0, 0, false, 1 ); // unknown 1 - y
bc.addCondition(boundary::east, condition_type::dirichlet, 0, 0, false, 2 ); // unknown 2 - z
// Surface forces
tmp << 0, 0, -90;
}

For the other test cases, please look at the file.

In the following snippet, the assembler is defined. Here, we can see that the gsExprAssembler is defined, typedefinitions for geometry map, variables, spaces and solutions are defined. Furthermore, the integration elements (quadrature points) are set according to the basis and an expression evaluator is constructed. Note that the gsExprAssembler::getSpace() function is called with the second argument being the dimension of he space. Furthermore, for the parameters and materials, variables are defined. Lastly, note the definition of the m2 object. This is a matrix with [1,1,2] on the diagonal because all the components (stresses, strains and their variations) are defined in Voight-notation.

gsExprAssembler<> A(1,1);
typedef gsExprAssembler<>::geometryMap geometryMap;
typedef gsExprAssembler<>::space space;
typedef gsExprAssembler<>::solution solution;
// Elements used for numerical integration
A.setIntegrationElements(dbasis);
gsExprEvaluator<> ev(A);
// Set the geometry map
geometryMap G = A.getMap(mp); // the last map counts
geometryMap defG = A.getMap(mp_def);
// Set the discretization space
space u = A.getSpace(dbasis, 3);
// Solution vector and solution variable
gsMatrix<> random;
solution u_sol = A.getSolution(u,random);
// gsFunctionExpr<> materialMat("1","0","0","0","1","0","0","0","1",3);
gsFunctionExpr<> E(util::to_string(E_modulus),3);
gsFunctionExpr<> nu(util::to_string(PoissonRatio),3);
gsMaterialMatrix<real_t> materialMat(mp, E, nu);
auto mm = A.getCoeff(materialMat); // evaluates in the parametric domain, but the class transforms E and nu to physical
gsFunctionExpr<> t(util::to_string(thickness), 3);
auto tt = A.getCoeff(t, G); // evaluates in the physical domain
gsFunctionExpr<> mult2t("1","0","0","0","1","0","0","0","2",3);
auto m2 = A.getCoeff(mult2t, G); // evaluates in the physical domain
gsConstantFunction<> force(tmp,3);
auto ff = A.getCoeff(force,G); // evaluates in the physical domain
gsSparseSolver<>::CGDiagonal solver;

Next up is the assembly of the system. First, all the necessary tensors (here, vectors by the Voight-notation) are defined. Dimensionality of the space is automatically taken into account. First, the weak Dirichlet and clamped boundary conditions are assembled according to the results of [5]. Then, the gsExprAssembler::assemble() function is called which assembles the linear system.

// Initialize the system
u.setup(bc, dirichlet::interpolation, 0);
A.initSystem();
gsInfo<<"Number of degrees of freedom: "<< A.numDofs() <<"\n"<<std::flush;
/*
We provide the following functions:
E_m membrane strain tensor.
E_m_der first variation of E_m.
E_m_der2 second variation of E_m MULTIPLIED BY S_m.
E_f flexural strain tensor.
E_f_der second variation of E_f.
E_f_der2 second variation of E_f MULTIPLIED BY S_f.
Where:
G the undeformed geometry,
defG the deformed geometry,
mm the material matrix,
m2 an auxillary matrix to multiply the last row of a tensor with 2
**/
// Membrane components
auto E_m = 0.5 * ( flat(jac(defG).tr()*jac(defG)) - flat(jac(G).tr()* jac(G)) ) ;
auto S_m = E_m * reshape(mm,3,3);
auto N = tt.val() * S_m;
auto E_m_der = flat( jac(defG).tr() * jac(u) ) ;
auto S_m_der = E_m_der * reshape(mm,3,3);
auto N_der = tt.val() * S_m_der;
auto E_m_der2 = flatdot( jac(u),jac(u).tr(), N );
// Flexural components
auto E_f = ( deriv2(G,sn(G).normalized().tr()) - deriv2(defG,sn(defG).normalized().tr()) ) * reshape(m2,3,3) ;
auto S_f = E_f * reshape(mm,3,3);
auto M = tt.val() * tt.val() * tt.val() / 12.0 * S_f;
auto E_f_der = -( deriv2(u,sn(defG).normalized().tr() ) + deriv2(defG,var1(u,defG) ) ) * reshape(m2,3,3);
auto S_f_der = E_f_der * reshape(mm,3,3);
auto M_der = tt.val() * tt.val() * tt.val() / 12.0 * S_f_der;
auto E_f_der2 = - (flatdot2( deriv2(u), var1(u,defG).tr(), M ).symmetrize() + var2(u,u,defG, M ));
auto F = ff;
// For Neumann (same for Dirichlet/Nitsche) conditions
auto g_N = A.getBdrFunction(G);
// auto g_N = ff;
real_t alpha_d = 1e3;
A.assembleBdr
(
bc.get("Weak Dirichlet")
,
alpha_d * u * u.tr()
,
alpha_d * (u * (defG - G) - u * (g_N) )
);
// For Neumann conditions
A.assembleBdr(bc.get("Neumann"), u * g_N * tv(G).norm() );
A.assemble(
(N_der * (E_m_der).tr() + M_der * (E_f_der).tr() ) * meas(G)
,
u * F * meas(G) + pressure * u * sn(defG).normalized() * meas(G)
);

The linear system can now be solved according tot the following. ALso, the deformed multipatch is constructed by adding the displacements - the solution to the system - to the original multipatch. Note that the solution vector contains the displacements of the control points of the original geometry. The solution object and the gsExprEvaluator are used for plotting. Also the stresses are evaluated at a point.

// solve system
solver.compute( A.matrix() );
gsMatrix<> solVector = solver.solve(A.rhs());
// update deformed patch
gsMatrix<> cc;
u_sol.setSolutionVector(solVector);
for ( size_t k =0; k!=mp_def.nPatches(); ++k) // Deform the geometry
{
// // extract deformed geometry
u_sol.extract(cc, k);
mp_def.patch(k).coefs() += cc; // defG points to mp_def, therefore updated
}
gsMatrix<> result;
u_sol.extractFull(result);

Now it is time for the nonlinear system assembly. The variational form depends on the deformed geometry which is the current geometry plus the computed displacements in the linear solve. In each iteration, the nonlinear variational form is assembled to a system of equations.

real_t residual = A.rhs().norm();
real_t residual0 = residual;
real_t residualOld = residual;
gsMatrix<> updateVector = solVector;
if (nonlinear)
{
index_t itMax = 25;
real_t tol = 1e-8;
for (index_t it = 0; it != itMax; ++it)
{
A.initSystem();
// assemble system
A.assemble(
(
N_der * E_m_der.tr()
+
E_m_der2
+
M_der * E_f_der.tr()
+
E_f_der2
) * meas(G)
-
pressure * u * var1(u,defG) .tr() * meas(G)
, u * F * meas(G) + pressure * u * sn(defG).normalized() * meas(G) - ( ( N * E_m_der.tr() + M * E_f_der.tr() ) * meas(G) ).tr()
);
// Neumann term
A.assembleBdr(bc.get("Neumann"), u * g_N * tv(G).norm() );
// Weak Dirichlet term
A.assembleBdr
(
bc.get("Weak Dirichlet")
,
alpha_d * u * u.tr()
,
-alpha_d * (u * (defG - G) - u * (g_N) )
);
// solve system
solver.compute( A.matrix() );
updateVector = solver.solve(A.rhs()); // this is the UPDATE
solVector += updateVector;
residual = A.rhs().norm();
gsInfo<<"Iteration: "<< it
<<", residue: "<< residual
<<", update norm: "<<updateVector.norm()
<<", log(Ri/R0): "<< math::log10(residualOld/residual0)
<<", log(Ri+1/R0): "<< math::log10(residual/residual0)
<<"\n";
residualOld = residual;
// update deformed patch
u_sol.setSolutionVector(updateVector);
for ( size_t k =0; k!=mp_def.nPatches(); ++k) // Deform the geometry
{
// // extract deformed geometry
u_sol.extract(cc, k);
mp_def.patch(k).coefs() += cc; // defG points to mp_def, therefore updated
}
if (residual < tol)
break;
}
}

When the solution to the nonlinear problem is known, the solution is constructed and the deformations are evaluated. The deformation is simply the difference between the deformed and the original geometries.

u_sol.setSolutionVector(solVector);
mp_def = mp;
for ( size_t k =0; k!=mp.nPatches(); ++k) // Deform the geometry
{
u_sol.extract(cc, k);
mp_def.patch(k).coefs() += cc; // defG points to mp_def, therefore updated
}
gsMultiPatch<> deformation = mp_def;
for (size_t k = 0; k != mp_def.nPatches(); ++k)
deformation.patch(k).coefs() -= mp.patch(k).coefs();
gsInfo <<"Maximum deformation coef: "
<< deformation.patch(0).coefs().colwise().maxCoeff() <<".\n";
gsInfo <<"Minimum deformation coef: "
<< deformation.patch(0).coefs().colwise().minCoeff() <<".\n";
gsInfo <<"Area (undeformed) = "<<ev.integral(meas(G))<<"\tArea (undeformed) = "<<ev.integral(meas(defG))<<"\n";

Lastly, the deformation can be evaluated on several points:

gsMatrix<> coords(2,1);
if (testCase==1)
coords<<0,0;
else if (testCase==2)
coords<<0.5,0;
else if (testCase==3)
coords<<0,0;
else if (testCase==4)
coords<<1,1;
else
coords<<0,0;
gsMatrix<> res(3,1);
mp.patch(0).eval_into(coords,res);
real_t x=res.at(0);
real_t y=res.at(1);
real_t z=res.at(2);
deformation.patch(0).eval_into(coords,res);
real_t ux=res.at(0);
real_t uy=res.at(1);
real_t uz=res.at(2);
gsInfo<<"Deformation on point ("<<x<<","<<y<<","<<z<<"):\n";
gsInfo<<std::setprecision(20)<<"("<<ux<<","<<uy<<","<<uz<<")"<<"\n";

And the solution can be plotted in Paraview:

if (plot)
{
gsField<> solField(mp, deformation);
gsInfo<<"Plotting in Paraview...\n";
gsWriteParaview<>( solField, "solution", 1000, true);
// gsFileManager::open("solution.pvd");
}

References

[1] J. Kiendl, K.-U. Bletzinger, J. Linhard, and R. Wüchner, “Isogeometric shell analysis with Kirchhoff–Love elements,” Comput. Methods Appl. Mech. Eng., vol. 198, no. 49–52, pp. 3902–3914, Nov. 2009.

[2] J. Kiendl, “Isogeometric analysis and shape optimal design of shell structures,” Technische Universität München, 2011.

[3] A. Goyal, “Isogeometric Shell Discretizations for Flexible Multibody Dynamics,” Technische Universität Kaiserslautern, 2015.

[4] H. M. Verhelst, “Modelling Wrinkling Behaviour of Large Floating Thin Offshore Structures: An application of Isogeometric Structural Analysis for Post-Buckling Analyses,” Delft University of Technology, 2019.

[5] A. J. Herrema, E. L. Johnson, D. Proserpio, M. C. H. Wu, J. Kiendl, and M.-C. Hsu, “Penalty coupling of non-matching isogeometric Kirchhoff–Love shell patches with application to composite wind turbine blades,” Comput. Methods Appl. Mech. Eng., vol. 346, pp. 810–840, Apr. 2019.

Contact

Author: Hugo Verhelst – h.m.v.nosp@m.erhe.nosp@m.lst@t.nosp@m.udel.nosp@m.ft.nl