xref: /libCEED/examples/solids/src/boundary.c (revision d612e6d93203685a17b8a19bf3945e10385cd221)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 /// @file
18 /// Boundary condition functions for solid mechanics example using PETSc
19 
20 #include "../elasticity.h"
21 
22 // -----------------------------------------------------------------------------
23 // Boundary Functions
24 // -----------------------------------------------------------------------------
25 // Note: If additional boundary conditions are added, an update is needed in
26 //         elasticity.h for the boundaryOptions variable.
27 
28 // BCMMS - boundary function
29 // Values on all points of the mesh is set based on given solution below
30 // for u[0], u[1], u[2]
31 PetscErrorCode BCMMS(PetscInt dim, PetscReal loadIncrement,
32                      const PetscReal coords[], PetscInt ncompu,
33                      PetscScalar *u, void *ctx) {
34   PetscScalar x = coords[0];
35   PetscScalar y = coords[1];
36   PetscScalar z = coords[2];
37 
38   PetscFunctionBeginUser;
39 
40   u[0] = exp(2*x)*sin(3*y)*cos(4*z) / 1e8 * loadIncrement;
41   u[1] = exp(3*y)*sin(4*z)*cos(2*x) / 1e8 * loadIncrement;
42   u[2] = exp(4*z)*sin(2*x)*cos(3*y) / 1e8 * loadIncrement;
43 
44   PetscFunctionReturn(0);
45 };
46 
47 // BCZero - fix boundary values at zero
48 PetscErrorCode BCZero(PetscInt dim, PetscReal loadIncrement,
49                       const PetscReal coords[], PetscInt ncompu,
50                       PetscScalar *u, void *ctx) {
51   PetscFunctionBeginUser;
52 
53   u[0] = 0;
54   u[1] = 0;
55   u[2] = 0;
56 
57   PetscFunctionReturn(0);
58 };
59 
60 // BCClampTranslate - translate boundary values at fraction of load increment
61 PetscErrorCode BCClampTranslate(PetscInt dim, PetscReal loadIncrement,
62                                 const PetscReal coords[], PetscInt ncompu,
63                                 PetscScalar *u, void *ctx) {
64   PetscScalar (*clampMax) = (PetscScalar(*))ctx;
65 
66   PetscFunctionBeginUser;
67 
68   u[0] = clampMax[0]*loadIncrement;
69   u[1] = clampMax[1]*loadIncrement;
70   u[2] = clampMax[2]*loadIncrement;
71 
72   PetscFunctionReturn(0);
73 };
74 
75 #ifndef M_PI
76 #  define M_PI    3.14159265358979323846
77 #endif
78 
79 // BCClampRotate - rotate boundary values at fraction of load increment
80 PetscErrorCode BCClampRotate(PetscInt dim, PetscReal loadIncrement,
81                              const PetscReal coords[], PetscInt ncompu,
82                              PetscScalar *u, void *ctx) {
83   PetscScalar x = coords[0];
84   PetscScalar y = coords[1];
85   PetscScalar z = coords[2];
86   PetscScalar (*clampMax) = (PetscScalar(*))ctx;
87 
88   PetscFunctionBeginUser;
89 
90   PetscScalar theta = clampMax[3]*M_PI*loadIncrement,
91               kx = clampMax[0], ky = clampMax[1], kz = clampMax[2];
92   PetscScalar c = cos(theta), s = sin(theta);
93 
94   u[0] = s*(-kz*y + ky*z) + (1-c)*(-ky*ky+kz*kz*x + kx*ky*y + kx*kz*z);
95   u[1] = s*(kz*x + -kx*z) + (1-c)*(kx*ky*x - (kx*kx+kz*kz)*y + ky*kz*z);
96   u[2] = s*(-ky*x + kx*y) + (1-c)*(kx*kz*x + ky*kz*y - (kx*kx+ky*ky)*z);
97 
98   PetscFunctionReturn(0);
99 };
100