xref: /libCEED/examples/solids/src/boundary.c (revision ea61e9ac44808524e4667c1525a05976f536c19c)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
33d8e8822SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53d8e8822SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d8e8822SJeremy L Thompson 
8ccaff030SJeremy L Thompson /// @file
9ccaff030SJeremy L Thompson /// Boundary condition functions for solid mechanics example using PETSc
10ccaff030SJeremy L Thompson 
115754ecacSJeremy L Thompson #include "../include/boundary.h"
12ccaff030SJeremy L Thompson 
13ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
14ccaff030SJeremy L Thompson // Boundary Functions
15ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
16*ea61e9acSJeremy L Thompson // Note: If additional boundary conditions are added, an update is needed in elasticity.h for the boundaryOptions variable.
17ccaff030SJeremy L Thompson 
18ccaff030SJeremy L Thompson // BCMMS - boundary function
19*ea61e9acSJeremy L Thompson // Values on all points of the mesh is set based on given solution below for u[0], u[1], u[2]
202b730f8bSJeremy L Thompson PetscErrorCode BCMMS(PetscInt dim, PetscReal load_increment, const PetscReal coords[], PetscInt num_comp_u, PetscScalar *u, void *ctx) {
21ccaff030SJeremy L Thompson   PetscScalar x = coords[0];
22ccaff030SJeremy L Thompson   PetscScalar y = coords[1];
23ccaff030SJeremy L Thompson   PetscScalar z = coords[2];
24ccaff030SJeremy L Thompson 
25ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
26ccaff030SJeremy L Thompson 
27d1d35e2fSjeremylt   u[0] = exp(2 * x) * sin(3 * y) * cos(4 * z) / 1e8 * load_increment;
28d1d35e2fSjeremylt   u[1] = exp(3 * y) * sin(4 * z) * cos(2 * x) / 1e8 * load_increment;
29d1d35e2fSjeremylt   u[2] = exp(4 * z) * sin(2 * x) * cos(3 * y) / 1e8 * load_increment;
30ccaff030SJeremy L Thompson 
31ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
32ccaff030SJeremy L Thompson };
33ccaff030SJeremy L Thompson 
3431dc5d86Sjeremylt #ifndef M_PI
3531dc5d86Sjeremylt #define M_PI 3.14159265358979323846
3631dc5d86Sjeremylt #endif
3731dc5d86Sjeremylt 
38*ea61e9acSJeremy L Thompson // BCClamp - fix boundary values with affine transformation at fraction of load increment
392b730f8bSJeremy L Thompson PetscErrorCode BCClamp(PetscInt dim, PetscReal load_increment, const PetscReal coords[], PetscInt num_comp_u, PetscScalar *u, void *ctx) {
4031dc5d86Sjeremylt   PetscScalar x          = coords[0];
4131dc5d86Sjeremylt   PetscScalar y          = coords[1];
4231dc5d86Sjeremylt   PetscScalar z          = coords[2];
4331dc5d86Sjeremylt   PetscScalar(*clampMax) = (PetscScalar(*))ctx;
4431dc5d86Sjeremylt 
4531dc5d86Sjeremylt   PetscFunctionBeginUser;
4672d03b64SArash Mehraban   PetscScalar
4772d03b64SArash Mehraban       // Translation vector
48d1d35e2fSjeremylt       lx = clampMax[0] * load_increment,
492b730f8bSJeremy L Thompson       ly = clampMax[1] * load_increment, lz = clampMax[2] * load_increment,
5072d03b64SArash Mehraban       // Normalized rotation axis
512b730f8bSJeremy L Thompson       kx = clampMax[3], ky = clampMax[4], kz = clampMax[5],
5272d03b64SArash Mehraban       // Rotation polynomial
532b730f8bSJeremy L Thompson       c_0 = clampMax[6] * M_PI, c_1 = clampMax[7] * M_PI, cx = kx * x + ky * y + kz * z,
5472d03b64SArash Mehraban       // Rotation magnitude
55d1d35e2fSjeremylt       theta     = (c_0 + c_1 * cx) * load_increment;
5631dc5d86Sjeremylt   PetscScalar c = cos(theta), s = sin(theta);
5731dc5d86Sjeremylt 
5856f0bea9Sjeremylt   u[0] = lx + s * (-kz * y + ky * z) + (1 - c) * (-(ky * ky + kz * kz) * x + kx * ky * y + kx * kz * z);
5956f0bea9Sjeremylt   u[1] = ly + s * (kz * x + -kx * z) + (1 - c) * (kx * ky * x + -(kx * kx + kz * kz) * y + ky * kz * z);
6056f0bea9Sjeremylt   u[2] = lz + s * (-ky * x + kx * y) + (1 - c) * (kx * kz * x + ky * kz * y + -(kx * kx + ky * ky) * z);
6131dc5d86Sjeremylt   PetscFunctionReturn(0);
6231dc5d86Sjeremylt };
63