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 #ifndef M_PI 48 # define M_PI 3.14159265358979323846 49 #endif 50 51 // BCClamp - fix boundary values with affine transformation at fraction of load 52 // increment 53 PetscErrorCode BCClamp(PetscInt dim, PetscReal loadIncrement, 54 const PetscReal coords[], PetscInt ncompu, 55 PetscScalar *u, void *ctx) { 56 PetscScalar x = coords[0]; 57 PetscScalar y = coords[1]; 58 PetscScalar z = coords[2]; 59 PetscScalar (*clampMax) = (PetscScalar(*))ctx; 60 61 PetscFunctionBeginUser; 62 63 PetscScalar lx = clampMax[0]*loadIncrement, ly = clampMax[1]*loadIncrement, 64 lz = clampMax[2]*loadIncrement, 65 theta = clampMax[6]*M_PI*loadIncrement, 66 kx = clampMax[3], ky = clampMax[4], kz = clampMax[5]; 67 PetscScalar c = cos(theta), s = sin(theta); 68 69 u[0] = lx + s*(-kz*y + ky*z) + (1-c)*(-(ky*ky+kz*kz)*x + kx*ky*y + kx*kz*z); 70 u[1] = ly + s*(kz*x + -kx*z) + (1-c)*(kx*ky*x + -(kx*kx+kz*kz)*y + ky*kz*z); 71 u[2] = lz + s*(-ky*x + kx*y) + (1-c)*(kx*kz*x + ky*kz*y + -(kx*kx+ky*ky)*z); 72 73 PetscFunctionReturn(0); 74 }; 75