xref: /libCEED/examples/solids/src/boundary.c (revision 00d548f6cce3ebb77e7917d79b50a2a368cad12e)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Boundary condition functions for solid mechanics example using PETSc
10 
11 #include "../include/boundary.h"
12 
13 // -----------------------------------------------------------------------------
14 // Boundary Functions
15 // -----------------------------------------------------------------------------
16 // Note: If additional boundary conditions are added, an update is needed in elasticity.h for the boundaryOptions variable.
17 
18 // BCMMS - boundary function
19 // Values on all points of the mesh is set based on given solution below for u[0], u[1], u[2]
20 PetscErrorCode BCMMS(PetscInt dim, PetscReal load_increment, const PetscReal coords[], PetscInt num_comp_u, PetscScalar *u, void *ctx) {
21   PetscScalar x = coords[0];
22   PetscScalar y = coords[1];
23   PetscScalar z = coords[2];
24 
25   PetscFunctionBeginUser;
26 
27   u[0] = exp(2 * x) * sin(3 * y) * cos(4 * z) / 1e8 * load_increment;
28   u[1] = exp(3 * y) * sin(4 * z) * cos(2 * x) / 1e8 * load_increment;
29   u[2] = exp(4 * z) * sin(2 * x) * cos(3 * y) / 1e8 * load_increment;
30 
31   PetscFunctionReturn(0);
32 };
33 
34 #ifndef M_PI
35 #define M_PI 3.14159265358979323846
36 #endif
37 
38 // BCClamp - fix boundary values with affine transformation at fraction of load increment
39 PetscErrorCode BCClamp(PetscInt dim, PetscReal load_increment, const PetscReal coords[], PetscInt num_comp_u, PetscScalar *u, void *ctx) {
40   PetscScalar x          = coords[0];
41   PetscScalar y          = coords[1];
42   PetscScalar z          = coords[2];
43   PetscScalar(*clampMax) = (PetscScalar(*))ctx;
44 
45   PetscFunctionBeginUser;
46   PetscScalar
47       // Translation vector
48       lx = clampMax[0] * load_increment,
49       ly = clampMax[1] * load_increment, lz = clampMax[2] * load_increment,
50       // Normalized rotation axis
51       kx = clampMax[3], ky = clampMax[4], kz = clampMax[5],
52       // Rotation polynomial
53       c_0 = clampMax[6] * M_PI, c_1 = clampMax[7] * M_PI, cx = kx * x + ky * y + kz * z,
54       // Rotation magnitude
55       theta     = (c_0 + c_1 * cx) * load_increment;
56   PetscScalar c = cos(theta), s = sin(theta);
57 
58   u[0] = lx + s * (-kz * y + ky * z) + (1 - c) * (-(ky * ky + kz * kz) * x + kx * ky * y + kx * kz * z);
59   u[1] = ly + s * (kz * x + -kx * z) + (1 - c) * (kx * ky * x + -(kx * kx + kz * kz) * y + ky * kz * z);
60   u[2] = lz + s * (-ky * x + kx * y) + (1 - c) * (kx * kz * x + ky * kz * y + -(kx * kx + ky * ky) * z);
61   PetscFunctionReturn(0);
62 };
63