xref: /libCEED/examples/nek/bps/bps.h (revision 4d537eea83dd2f1011134892241345b6ac537f56)
1*4d537eeaSYohann // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC.
2*4d537eeaSYohann // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707.
3*4d537eeaSYohann // All Rights reserved. See files LICENSE and NOTICE for details.
4*4d537eeaSYohann //
5*4d537eeaSYohann // This file is part of CEED, a collection of benchmarks, miniapps, software
6*4d537eeaSYohann // libraries and APIs for efficient high-order finite element and spectral
7*4d537eeaSYohann // element discretizations for exascale applications. For more information and
8*4d537eeaSYohann // source code availability see http://github.com/ceed.
9*4d537eeaSYohann //
10*4d537eeaSYohann // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11*4d537eeaSYohann // a collaborative effort of two U.S. Department of Energy organizations (Office
12*4d537eeaSYohann // of Science and the National Nuclear Security Administration) responsible for
13*4d537eeaSYohann // the planning and preparation of a capable exascale ecosystem, including
14*4d537eeaSYohann // software, applications, hardware, advanced system engineering and early
15*4d537eeaSYohann // testbed platforms, in support of the nation's exascale computing imperative.
16*4d537eeaSYohann 
17*4d537eeaSYohann #ifndef __CUDACC__
18*4d537eeaSYohann #  include <math.h>
19*4d537eeaSYohann #endif
20*4d537eeaSYohann 
21*4d537eeaSYohann // *****************************************************************************
22*4d537eeaSYohann //   BP 1
23*4d537eeaSYohann // *****************************************************************************
24*4d537eeaSYohann CEED_QFUNCTION(masssetupf)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
25*4d537eeaSYohann   CeedScalar *rho = out[0], *rhs = out[1];
26*4d537eeaSYohann   const CeedScalar *x = in[0];
27*4d537eeaSYohann   const CeedScalar *J = in[1];
28*4d537eeaSYohann   const CeedScalar *w = in[2];
29*4d537eeaSYohann   for (CeedInt i=0; i<Q; i++) {
30*4d537eeaSYohann     CeedScalar det = (J[i+Q*0]*(J[i+Q*4]*J[i+Q*8] - J[i+Q*5]*J[i+Q*7]) -
31*4d537eeaSYohann                       J[i+Q*1]*(J[i+Q*3]*J[i+Q*8] - J[i+Q*5]*J[i+Q*6]) +
32*4d537eeaSYohann                       J[i+Q*2]*(J[i+Q*3]*J[i+Q*7] - J[i+Q*4]*J[i+Q*6]));
33*4d537eeaSYohann     rho[i] = det * w[i];
34*4d537eeaSYohann     rhs[i] = rho[i] * w[i] *
35*4d537eeaSYohann                sqrt(x[i]*x[i] + x[i+Q]*x[i+Q] + x[i+2*Q]*x[i+2*Q]);
36*4d537eeaSYohann   }
37*4d537eeaSYohann   return 0;
38*4d537eeaSYohann }
39*4d537eeaSYohann 
40*4d537eeaSYohann CEED_QFUNCTION int massf(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
41*4d537eeaSYohann   const CeedScalar *u = in[0];
42*4d537eeaSYohann   const CeedScalar *rho = in[1];
43*4d537eeaSYohann   CeedScalar *v = out[0];
44*4d537eeaSYohann   for (CeedInt i=0; i<Q; i++) {
45*4d537eeaSYohann     v[i] = rho[i] * u[i];
46*4d537eeaSYohann   }
47*4d537eeaSYohann   return 0;
48*4d537eeaSYohann }
49*4d537eeaSYohann // *****************************************************************************
50*4d537eeaSYohann //   BP 3
51*4d537eeaSYohann // *****************************************************************************
52*4d537eeaSYohann CEED_QFUNCTION(diffsetupf)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
53*4d537eeaSYohann   #ifndef M_PI
54*4d537eeaSYohann   #define M_PI    3.14159265358979323846
55*4d537eeaSYohann   #endif
56*4d537eeaSYohann   const CeedScalar *x = in[0];
57*4d537eeaSYohann   const CeedScalar *J = in[1];
58*4d537eeaSYohann   const CeedScalar *w = in[2];
59*4d537eeaSYohann   CeedScalar *qd = out[0], *rhs = out[1];
60*4d537eeaSYohann   for (CeedInt i=0; i<Q; i++) {
61*4d537eeaSYohann     const CeedScalar J11 = J[i+Q*0];
62*4d537eeaSYohann     const CeedScalar J21 = J[i+Q*1];
63*4d537eeaSYohann     const CeedScalar J31 = J[i+Q*2];
64*4d537eeaSYohann     const CeedScalar J12 = J[i+Q*3];
65*4d537eeaSYohann     const CeedScalar J22 = J[i+Q*4];
66*4d537eeaSYohann     const CeedScalar J32 = J[i+Q*5];
67*4d537eeaSYohann     const CeedScalar J13 = J[i+Q*6];
68*4d537eeaSYohann     const CeedScalar J23 = J[i+Q*7];
69*4d537eeaSYohann     const CeedScalar J33 = J[i+Q*8];
70*4d537eeaSYohann     const CeedScalar A11 = J22*J33 - J23*J32;
71*4d537eeaSYohann     const CeedScalar A12 = J13*J32 - J12*J33;
72*4d537eeaSYohann     const CeedScalar A13 = J12*J23 - J13*J22;
73*4d537eeaSYohann     const CeedScalar A21 = J23*J31 - J21*J33;
74*4d537eeaSYohann     const CeedScalar A22 = J11*J33 - J13*J31;
75*4d537eeaSYohann     const CeedScalar A23 = J13*J21 - J11*J23;
76*4d537eeaSYohann     const CeedScalar A31 = J21*J32 - J22*J31;
77*4d537eeaSYohann     const CeedScalar A32 = J12*J31 - J11*J32;
78*4d537eeaSYohann     const CeedScalar A33 = J11*J22 - J12*J21;
79*4d537eeaSYohann     const CeedScalar qw = w[i] / (J11*A11 + J21*A12 + J31*A13);
80*4d537eeaSYohann     qd[i+Q*0] = qw * (A11*A11 + A12*A12 + A13*A13);
81*4d537eeaSYohann     qd[i+Q*1] = qw * (A11*A21 + A12*A22 + A13*A23);
82*4d537eeaSYohann     qd[i+Q*2] = qw * (A11*A31 + A12*A32 + A13*A33);
83*4d537eeaSYohann     qd[i+Q*3] = qw * (A21*A21 + A22*A22 + A23*A23);
84*4d537eeaSYohann     qd[i+Q*4] = qw * (A21*A31 + A22*A32 + A23*A33);
85*4d537eeaSYohann     qd[i+Q*5] = qw * (A31*A31 + A32*A32 + A33*A33);
86*4d537eeaSYohann     const CeedScalar c[3] = { 0, 1., 2. };
87*4d537eeaSYohann     const CeedScalar k[3] = { 1., 2., 3. };
88*4d537eeaSYohann     const CeedScalar rho = w[i] * (J11*A11 + J21*A12 + J31*A13);
89*4d537eeaSYohann     rhs[i] = rho * M_PI*M_PI * (k[0]*k[0] + k[1]*k[1] + k[2]*k[2]) *
90*4d537eeaSYohann                sin(M_PI*(c[0] + k[0]*x[i+Q*0])) *
91*4d537eeaSYohann                sin(M_PI*(c[1] + k[1]*x[i+Q*1])) *
92*4d537eeaSYohann                sin(M_PI*(c[2] + k[2]*x[i+Q*2]));
93*4d537eeaSYohann   }
94*4d537eeaSYohann   return 0;
95*4d537eeaSYohann }
96*4d537eeaSYohann 
97*4d537eeaSYohann CEED_QFUNCTION int diffusionf(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
98*4d537eeaSYohann   const CeedScalar *ug = in[0];
99*4d537eeaSYohann   const CeedScalar *qd = in[1];
100*4d537eeaSYohann   CeedScalar *vg = out[0];
101*4d537eeaSYohann   for (CeedInt i=0; i<Q; i++) {
102*4d537eeaSYohann     const CeedScalar ug0 = ug[i+Q*0];
103*4d537eeaSYohann     const CeedScalar ug1 = ug[i+Q*1];
104*4d537eeaSYohann     const CeedScalar ug2 = ug[i+Q*2];
105*4d537eeaSYohann     vg[i+Q*0] = qd[i+Q*0]*ug0 + qd[i+Q*1]*ug1 + qd[i+Q*2]*ug2;
106*4d537eeaSYohann     vg[i+Q*1] = qd[i+Q*1]*ug0 + qd[i+Q*3]*ug1 + qd[i+Q*4]*ug2;
107*4d537eeaSYohann     vg[i+Q*2] = qd[i+Q*2]*ug0 + qd[i+Q*4]*ug1 + qd[i+Q*5]*ug2;
108*4d537eeaSYohann   }
109*4d537eeaSYohann   return 0;
110*4d537eeaSYohann }
111