1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscsys.h>
1549aac155SJeremy L Thompson
16ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
17ccaff030SJeremy L Thompson // Boundary Functions
18ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
19ea61e9acSJeremy L Thompson // Note: If additional boundary conditions are added, an update is needed in elasticity.h for the boundaryOptions variable.
20ccaff030SJeremy L Thompson
21ccaff030SJeremy L Thompson // BCMMS - boundary function
22ea61e9acSJeremy L Thompson // Values on all points of the mesh is set based on given solution below for u[0], u[1], u[2]
BCMMS(PetscInt dim,PetscReal load_increment,const PetscReal coords[],PetscInt num_comp_u,PetscScalar * u,void * ctx)232b730f8bSJeremy L Thompson PetscErrorCode BCMMS(PetscInt dim, PetscReal load_increment, const PetscReal coords[], PetscInt num_comp_u, PetscScalar *u, void *ctx) {
24ccaff030SJeremy L Thompson PetscScalar x = coords[0];
25ccaff030SJeremy L Thompson PetscScalar y = coords[1];
26ccaff030SJeremy L Thompson PetscScalar z = coords[2];
27ccaff030SJeremy L Thompson
28ccaff030SJeremy L Thompson PetscFunctionBeginUser;
29ccaff030SJeremy L Thompson
30d1d35e2fSjeremylt u[0] = exp(2 * x) * sin(3 * y) * cos(4 * z) / 1e8 * load_increment;
31d1d35e2fSjeremylt u[1] = exp(3 * y) * sin(4 * z) * cos(2 * x) / 1e8 * load_increment;
32d1d35e2fSjeremylt u[2] = exp(4 * z) * sin(2 * x) * cos(3 * y) / 1e8 * load_increment;
33ccaff030SJeremy L Thompson
34ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
35ccaff030SJeremy L Thompson };
36ccaff030SJeremy L Thompson
3731dc5d86Sjeremylt #ifndef M_PI
3831dc5d86Sjeremylt #define M_PI 3.14159265358979323846
3931dc5d86Sjeremylt #endif
4031dc5d86Sjeremylt
41ea61e9acSJeremy L Thompson // BCClamp - fix boundary values with affine transformation at fraction of load increment
BCClamp(PetscInt dim,PetscReal load_increment,const PetscReal coords[],PetscInt num_comp_u,PetscScalar * u,void * ctx)422b730f8bSJeremy L Thompson PetscErrorCode BCClamp(PetscInt dim, PetscReal load_increment, const PetscReal coords[], PetscInt num_comp_u, PetscScalar *u, void *ctx) {
4331dc5d86Sjeremylt PetscScalar x = coords[0];
4431dc5d86Sjeremylt PetscScalar y = coords[1];
4531dc5d86Sjeremylt PetscScalar z = coords[2];
4631dc5d86Sjeremylt PetscScalar(*clampMax) = (PetscScalar(*))ctx;
4731dc5d86Sjeremylt
4831dc5d86Sjeremylt PetscFunctionBeginUser;
4972d03b64SArash Mehraban PetscScalar
5072d03b64SArash Mehraban // Translation vector
51d1d35e2fSjeremylt lx = clampMax[0] * load_increment,
522b730f8bSJeremy L Thompson ly = clampMax[1] * load_increment, lz = clampMax[2] * load_increment,
5372d03b64SArash Mehraban // Normalized rotation axis
542b730f8bSJeremy L Thompson kx = clampMax[3], ky = clampMax[4], kz = clampMax[5],
5572d03b64SArash Mehraban // Rotation polynomial
562b730f8bSJeremy L Thompson c_0 = clampMax[6] * M_PI, c_1 = clampMax[7] * M_PI, cx = kx * x + ky * y + kz * z,
5772d03b64SArash Mehraban // Rotation magnitude
58d1d35e2fSjeremylt theta = (c_0 + c_1 * cx) * load_increment;
5931dc5d86Sjeremylt PetscScalar c = cos(theta), s = sin(theta);
6031dc5d86Sjeremylt
6156f0bea9Sjeremylt u[0] = lx + s * (-kz * y + ky * z) + (1 - c) * (-(ky * ky + kz * kz) * x + kx * ky * y + kx * kz * z);
6256f0bea9Sjeremylt u[1] = ly + s * (kz * x + -kx * z) + (1 - c) * (kx * ky * x + -(kx * kx + kz * kz) * y + ky * kz * z);
6356f0bea9Sjeremylt u[2] = lz + s * (-ky * x + kx * y) + (1 - c) * (kx * kz * x + ky * kz * y + -(kx * kx + ky * ky) * z);
64ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
6531dc5d86Sjeremylt };
66