xref: /libCEED/examples/solids/src/boundary.c (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*3d8e8822SJeremy L Thompson //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5*3d8e8822SJeremy L Thompson //
6*3d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7*3d8e8822SJeremy 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 // -----------------------------------------------------------------------------
16ccaff030SJeremy L Thompson // Note: If additional boundary conditions are added, an update is needed in
17ccaff030SJeremy L Thompson //         elasticity.h for the boundaryOptions variable.
18ccaff030SJeremy L Thompson 
19ccaff030SJeremy L Thompson // BCMMS - boundary function
20ccaff030SJeremy L Thompson // Values on all points of the mesh is set based on given solution below
21ccaff030SJeremy L Thompson // for u[0], u[1], u[2]
22d1d35e2fSjeremylt PetscErrorCode BCMMS(PetscInt dim, PetscReal load_increment,
23d1d35e2fSjeremylt                      const PetscReal coords[], PetscInt num_comp_u,
24ccaff030SJeremy L Thompson                      PetscScalar *u, void *ctx) {
25ccaff030SJeremy L Thompson   PetscScalar x = coords[0];
26ccaff030SJeremy L Thompson   PetscScalar y = coords[1];
27ccaff030SJeremy L Thompson   PetscScalar z = coords[2];
28ccaff030SJeremy L Thompson 
29ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
30ccaff030SJeremy L Thompson 
31d1d35e2fSjeremylt   u[0] = exp(2*x)*sin(3*y)*cos(4*z) / 1e8 * load_increment;
32d1d35e2fSjeremylt   u[1] = exp(3*y)*sin(4*z)*cos(2*x) / 1e8 * load_increment;
33d1d35e2fSjeremylt   u[2] = exp(4*z)*sin(2*x)*cos(3*y) / 1e8 * load_increment;
34ccaff030SJeremy L Thompson 
35ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
36ccaff030SJeremy L Thompson };
37ccaff030SJeremy L Thompson 
3831dc5d86Sjeremylt #ifndef M_PI
3931dc5d86Sjeremylt #  define M_PI    3.14159265358979323846
4031dc5d86Sjeremylt #endif
4131dc5d86Sjeremylt 
42d642641fSjeremylt // BCClamp - fix boundary values with affine transformation at fraction of load
43d642641fSjeremylt //   increment
44d1d35e2fSjeremylt PetscErrorCode BCClamp(PetscInt dim, PetscReal load_increment,
45d1d35e2fSjeremylt                        const PetscReal coords[], PetscInt num_comp_u,
4631dc5d86Sjeremylt                        PetscScalar *u, void *ctx) {
4731dc5d86Sjeremylt   PetscScalar x = coords[0];
4831dc5d86Sjeremylt   PetscScalar y = coords[1];
4931dc5d86Sjeremylt   PetscScalar z = coords[2];
5031dc5d86Sjeremylt   PetscScalar (*clampMax) = (PetscScalar(*))ctx;
5131dc5d86Sjeremylt 
5231dc5d86Sjeremylt   PetscFunctionBeginUser;
5372d03b64SArash Mehraban   PetscScalar
5472d03b64SArash Mehraban   // Translation vector
55d1d35e2fSjeremylt   lx = clampMax[0]*load_increment,
56d1d35e2fSjeremylt   ly = clampMax[1]*load_increment,
57d1d35e2fSjeremylt   lz = clampMax[2]*load_increment,
5872d03b64SArash Mehraban   // Normalized rotation axis
5972d03b64SArash Mehraban   kx = clampMax[3],
6072d03b64SArash Mehraban   ky = clampMax[4],
6172d03b64SArash Mehraban   kz = clampMax[5],
6272d03b64SArash Mehraban   // Rotation polynomial
6372d03b64SArash Mehraban   c_0 = clampMax[6] * M_PI,
6472d03b64SArash Mehraban   c_1 = clampMax[7] * M_PI,
6572d03b64SArash Mehraban   cx = kx * x + ky * y + kz * z,
6672d03b64SArash Mehraban   // Rotation magnitude
67d1d35e2fSjeremylt   theta = (c_0 + c_1 * cx) * load_increment;
6831dc5d86Sjeremylt   PetscScalar c = cos(theta), s = sin(theta);
6931dc5d86Sjeremylt 
7056f0bea9Sjeremylt   u[0] = lx + s*(-kz*y + ky*z) + (1-c)*(-(ky*ky+kz*kz)*x + kx*ky*y + kx*kz*z);
7156f0bea9Sjeremylt   u[1] = ly + s*(kz*x + -kx*z) + (1-c)*(kx*ky*x + -(kx*kx+kz*kz)*y + ky*kz*z);
7256f0bea9Sjeremylt   u[2] = lz + s*(-ky*x + kx*y) + (1-c)*(kx*kz*x + ky*kz*y + -(kx*kx+ky*ky)*z);
7331dc5d86Sjeremylt   PetscFunctionReturn(0);
7431dc5d86Sjeremylt };
75