xref: /libCEED/examples/mfem/bp3.h (revision a57ca787a7f60d7c960eb3632298c2a0ffab656b)
14d537eeaSYohann // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC.
24d537eeaSYohann // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707.
34d537eeaSYohann // All Rights reserved. See files LICENSE and NOTICE for details.
44d537eeaSYohann //
54d537eeaSYohann // This file is part of CEED, a collection of benchmarks, miniapps, software
64d537eeaSYohann // libraries and APIs for efficient high-order finite element and spectral
74d537eeaSYohann // element discretizations for exascale applications. For more information and
84d537eeaSYohann // source code availability see http://github.com/ceed.
94d537eeaSYohann //
104d537eeaSYohann // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
114d537eeaSYohann // a collaborative effort of two U.S. Department of Energy organizations (Office
124d537eeaSYohann // of Science and the National Nuclear Security Administration) responsible for
134d537eeaSYohann // the planning and preparation of a capable exascale ecosystem, including
144d537eeaSYohann // software, applications, hardware, advanced system engineering and early
154d537eeaSYohann // testbed platforms, in support of the nation's exascale computing imperative.
164d537eeaSYohann 
17*a57ca787Svaleriabarra #ifndef bp3_h
18*a57ca787Svaleriabarra #define bp3_h
19*a57ca787Svaleriabarra #include <ceed.h>
204d537eeaSYohann 
214d537eeaSYohann /// A structure used to pass additional data to f_build_diff and f_apply_diff
224d537eeaSYohann struct BuildContext { CeedInt dim, space_dim; };
234d537eeaSYohann 
244d537eeaSYohann /// libCEED Q-function for building quadrature data for a diffusion operator
254d537eeaSYohann CEED_QFUNCTION(f_build_diff)(void *ctx, const CeedInt Q,
264d537eeaSYohann                              const CeedScalar *const *in, CeedScalar *const *out) {
274d537eeaSYohann   BuildContext *bc = (BuildContext *)ctx;
284d537eeaSYohann   // in[0] is Jacobians with shape [dim, nc=dim, Q]
294d537eeaSYohann   // in[1] is quadrature weights, size (Q)
304d537eeaSYohann   //
31a2fa7910SValeria Barra   // At every quadrature point, compute w/det(J).adj(J).adj(J)^T and store
324d537eeaSYohann   // the symmetric part of the result.
33a2fa7910SValeria Barra   const CeedScalar *J = in[0], *w = in[1];
34a2fa7910SValeria Barra   CeedScalar *qdata = out[0];
35ee07ded2SValeria Barra 
364d537eeaSYohann   switch (bc->dim + 10*bc->space_dim) {
374d537eeaSYohann   case 11:
38ee07ded2SValeria Barra     // Quadrature Point Loop
39ee07ded2SValeria Barra     CeedPragmaSIMD
404d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
41a2fa7910SValeria Barra       qdata[i] = w[i] / J[i];
424d537eeaSYohann     }
434d537eeaSYohann     break;
444d537eeaSYohann   case 22:
45ee07ded2SValeria Barra     // Quadrature Point Loop
46ee07ded2SValeria Barra     CeedPragmaSIMD
474d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
48a2fa7910SValeria Barra       // J: 0 2   qdata: 0 2   adj(J):  J22 -J12
49288c0443SJeremy L Thompson       //    1 3          2 1           -J21  J11
504d537eeaSYohann       const CeedScalar J11 = J[i+Q*0];
514d537eeaSYohann       const CeedScalar J21 = J[i+Q*1];
524d537eeaSYohann       const CeedScalar J12 = J[i+Q*2];
534d537eeaSYohann       const CeedScalar J22 = J[i+Q*3];
54a2fa7910SValeria Barra       const CeedScalar qw = w[i] / (J11*J22 - J21*J12);
55a2fa7910SValeria Barra       qdata[i+Q*0] =   qw * (J12*J12 + J22*J22);
56a2fa7910SValeria Barra       qdata[i+Q*1] =   qw * (J11*J11 + J21*J21);
57a2fa7910SValeria Barra       qdata[i+Q*2] = - qw * (J11*J12 + J21*J22);
584d537eeaSYohann     }
594d537eeaSYohann     break;
604d537eeaSYohann   case 33:
61ee07ded2SValeria Barra     // Quadrature Point Loop
62ee07ded2SValeria Barra     CeedPragmaSIMD
634d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
64a2fa7910SValeria Barra       // J: 0 3 6   qdata: 0 5 4
65288c0443SJeremy L Thompson       //    1 4 7          5 1 3
66288c0443SJeremy L Thompson       //    2 5 8          4 3 2
674d537eeaSYohann       const CeedScalar J11 = J[i+Q*0];
684d537eeaSYohann       const CeedScalar J21 = J[i+Q*1];
694d537eeaSYohann       const CeedScalar J31 = J[i+Q*2];
704d537eeaSYohann       const CeedScalar J12 = J[i+Q*3];
714d537eeaSYohann       const CeedScalar J22 = J[i+Q*4];
724d537eeaSYohann       const CeedScalar J32 = J[i+Q*5];
734d537eeaSYohann       const CeedScalar J13 = J[i+Q*6];
744d537eeaSYohann       const CeedScalar J23 = J[i+Q*7];
754d537eeaSYohann       const CeedScalar J33 = J[i+Q*8];
764d537eeaSYohann       const CeedScalar A11 = J22*J33 - J23*J32;
774d537eeaSYohann       const CeedScalar A12 = J13*J32 - J12*J33;
784d537eeaSYohann       const CeedScalar A13 = J12*J23 - J13*J22;
794d537eeaSYohann       const CeedScalar A21 = J23*J31 - J21*J33;
804d537eeaSYohann       const CeedScalar A22 = J11*J33 - J13*J31;
814d537eeaSYohann       const CeedScalar A23 = J13*J21 - J11*J23;
824d537eeaSYohann       const CeedScalar A31 = J21*J32 - J22*J31;
834d537eeaSYohann       const CeedScalar A32 = J12*J31 - J11*J32;
844d537eeaSYohann       const CeedScalar A33 = J11*J22 - J12*J21;
85a2fa7910SValeria Barra       const CeedScalar qw = w[i] / (J11*A11 + J21*A12 + J31*A13);
86a2fa7910SValeria Barra       qdata[i+Q*0] = qw * (A11*A11 + A12*A12 + A13*A13);
87a2fa7910SValeria Barra       qdata[i+Q*1] = qw * (A21*A21 + A22*A22 + A23*A23);
88a2fa7910SValeria Barra       qdata[i+Q*2] = qw * (A31*A31 + A32*A32 + A33*A33);
89a2fa7910SValeria Barra       qdata[i+Q*3] = qw * (A21*A31 + A22*A32 + A23*A33);
90a2fa7910SValeria Barra       qdata[i+Q*4] = qw * (A11*A31 + A12*A32 + A13*A33);
91a2fa7910SValeria Barra       qdata[i+Q*5] = qw * (A11*A21 + A12*A22 + A13*A23);
924d537eeaSYohann     }
934d537eeaSYohann     break;
944d537eeaSYohann   }
954d537eeaSYohann   return 0;
964d537eeaSYohann }
974d537eeaSYohann 
984d537eeaSYohann /// libCEED Q-function for applying a diff operator
994d537eeaSYohann CEED_QFUNCTION(f_apply_diff)(void *ctx, const CeedInt Q,
1004d537eeaSYohann                              const CeedScalar *const *in, CeedScalar *const *out) {
1014d537eeaSYohann   BuildContext *bc = (BuildContext *)ctx;
1024d537eeaSYohann   // in[0], out[0] have shape [dim, nc=1, Q]
103a2fa7910SValeria Barra   const CeedScalar *ug = in[0], *qdata = in[1];
1044d537eeaSYohann   CeedScalar *vg = out[0];
105ee07ded2SValeria Barra 
1064d537eeaSYohann   switch (bc->dim) {
1074d537eeaSYohann   case 1:
108ee07ded2SValeria Barra     // Quadrature Point Loop
109ee07ded2SValeria Barra     CeedPragmaSIMD
1104d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
111a2fa7910SValeria Barra       vg[i] = ug[i] * qdata[i];
1124d537eeaSYohann     }
1134d537eeaSYohann     break;
1144d537eeaSYohann   case 2:
115ee07ded2SValeria Barra     // Quadrature Point Loop
116ee07ded2SValeria Barra     CeedPragmaSIMD
1174d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
1184d537eeaSYohann       const CeedScalar ug0 = ug[i+Q*0];
1194d537eeaSYohann       const CeedScalar ug1 = ug[i+Q*1];
120a2fa7910SValeria Barra       vg[i+Q*0] = qdata[i+Q*0]*ug0 + qdata[i+Q*2]*ug1;
121a2fa7910SValeria Barra       vg[i+Q*1] = qdata[i+Q*2]*ug0 + qdata[i+Q*1]*ug1;
1224d537eeaSYohann     }
1234d537eeaSYohann     break;
1244d537eeaSYohann   case 3:
125ee07ded2SValeria Barra     // Quadrature Point Loop
126ee07ded2SValeria Barra     CeedPragmaSIMD
1274d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
1284d537eeaSYohann       const CeedScalar ug0 = ug[i+Q*0];
1294d537eeaSYohann       const CeedScalar ug1 = ug[i+Q*1];
1304d537eeaSYohann       const CeedScalar ug2 = ug[i+Q*2];
131a2fa7910SValeria Barra       vg[i+Q*0] = qdata[i+Q*0]*ug0 + qdata[i+Q*5]*ug1 + qdata[i+Q*4]*ug2;
132a2fa7910SValeria Barra       vg[i+Q*1] = qdata[i+Q*5]*ug0 + qdata[i+Q*1]*ug1 + qdata[i+Q*3]*ug2;
133a2fa7910SValeria Barra       vg[i+Q*2] = qdata[i+Q*4]*ug0 + qdata[i+Q*3]*ug1 + qdata[i+Q*2]*ug2;
1344d537eeaSYohann     }
1354d537eeaSYohann     break;
1364d537eeaSYohann   }
1374d537eeaSYohann   return 0;
1384d537eeaSYohann }
139288c0443SJeremy L Thompson 
140*a57ca787Svaleriabarra #endif // bp3_h
141