xref: /libCEED/examples/mfem/bp3.h (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.
34d537eeaSYohann //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
54d537eeaSYohann //
6*3d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
74d537eeaSYohann 
8a57ca787Svaleriabarra #ifndef bp3_h
9a57ca787Svaleriabarra #define bp3_h
104d537eeaSYohann 
114d537eeaSYohann /// A structure used to pass additional data to f_build_diff and f_apply_diff
124d537eeaSYohann struct BuildContext { CeedInt dim, space_dim; };
134d537eeaSYohann 
144d537eeaSYohann /// libCEED Q-function for building quadrature data for a diffusion operator
154d537eeaSYohann CEED_QFUNCTION(f_build_diff)(void *ctx, const CeedInt Q,
164d537eeaSYohann                              const CeedScalar *const *in, CeedScalar *const *out) {
174d537eeaSYohann   BuildContext *bc = (BuildContext *)ctx;
184d537eeaSYohann   // in[0] is Jacobians with shape [dim, nc=dim, Q]
194d537eeaSYohann   // in[1] is quadrature weights, size (Q)
204d537eeaSYohann   //
21a2fa7910SValeria Barra   // At every quadrature point, compute w/det(J).adj(J).adj(J)^T and store
224d537eeaSYohann   // the symmetric part of the result.
23a2fa7910SValeria Barra   const CeedScalar *J = in[0], *w = in[1];
24a2fa7910SValeria Barra   CeedScalar *qdata = out[0];
25ee07ded2SValeria Barra 
264d537eeaSYohann   switch (bc->dim + 10*bc->space_dim) {
274d537eeaSYohann   case 11:
28ee07ded2SValeria Barra     // Quadrature Point Loop
29ee07ded2SValeria Barra     CeedPragmaSIMD
304d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
31a2fa7910SValeria Barra       qdata[i] = w[i] / J[i];
324d537eeaSYohann     }
334d537eeaSYohann     break;
344d537eeaSYohann   case 22:
35ee07ded2SValeria Barra     // Quadrature Point Loop
36ee07ded2SValeria Barra     CeedPragmaSIMD
374d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
38a2fa7910SValeria Barra       // J: 0 2   qdata: 0 2   adj(J):  J22 -J12
39288c0443SJeremy L Thompson       //    1 3          2 1           -J21  J11
404d537eeaSYohann       const CeedScalar J11 = J[i+Q*0];
414d537eeaSYohann       const CeedScalar J21 = J[i+Q*1];
424d537eeaSYohann       const CeedScalar J12 = J[i+Q*2];
434d537eeaSYohann       const CeedScalar J22 = J[i+Q*3];
44a2fa7910SValeria Barra       const CeedScalar qw = w[i] / (J11*J22 - J21*J12);
45a2fa7910SValeria Barra       qdata[i+Q*0] =   qw * (J12*J12 + J22*J22);
46a2fa7910SValeria Barra       qdata[i+Q*1] =   qw * (J11*J11 + J21*J21);
47a2fa7910SValeria Barra       qdata[i+Q*2] = - qw * (J11*J12 + J21*J22);
484d537eeaSYohann     }
494d537eeaSYohann     break;
504d537eeaSYohann   case 33:
51ee07ded2SValeria Barra     // Quadrature Point Loop
52ee07ded2SValeria Barra     CeedPragmaSIMD
534d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
54a2fa7910SValeria Barra       // J: 0 3 6   qdata: 0 5 4
55288c0443SJeremy L Thompson       //    1 4 7          5 1 3
56288c0443SJeremy L Thompson       //    2 5 8          4 3 2
574d537eeaSYohann       const CeedScalar J11 = J[i+Q*0];
584d537eeaSYohann       const CeedScalar J21 = J[i+Q*1];
594d537eeaSYohann       const CeedScalar J31 = J[i+Q*2];
604d537eeaSYohann       const CeedScalar J12 = J[i+Q*3];
614d537eeaSYohann       const CeedScalar J22 = J[i+Q*4];
624d537eeaSYohann       const CeedScalar J32 = J[i+Q*5];
634d537eeaSYohann       const CeedScalar J13 = J[i+Q*6];
644d537eeaSYohann       const CeedScalar J23 = J[i+Q*7];
654d537eeaSYohann       const CeedScalar J33 = J[i+Q*8];
664d537eeaSYohann       const CeedScalar A11 = J22*J33 - J23*J32;
674d537eeaSYohann       const CeedScalar A12 = J13*J32 - J12*J33;
684d537eeaSYohann       const CeedScalar A13 = J12*J23 - J13*J22;
694d537eeaSYohann       const CeedScalar A21 = J23*J31 - J21*J33;
704d537eeaSYohann       const CeedScalar A22 = J11*J33 - J13*J31;
714d537eeaSYohann       const CeedScalar A23 = J13*J21 - J11*J23;
724d537eeaSYohann       const CeedScalar A31 = J21*J32 - J22*J31;
734d537eeaSYohann       const CeedScalar A32 = J12*J31 - J11*J32;
744d537eeaSYohann       const CeedScalar A33 = J11*J22 - J12*J21;
75a2fa7910SValeria Barra       const CeedScalar qw = w[i] / (J11*A11 + J21*A12 + J31*A13);
76a2fa7910SValeria Barra       qdata[i+Q*0] = qw * (A11*A11 + A12*A12 + A13*A13);
77a2fa7910SValeria Barra       qdata[i+Q*1] = qw * (A21*A21 + A22*A22 + A23*A23);
78a2fa7910SValeria Barra       qdata[i+Q*2] = qw * (A31*A31 + A32*A32 + A33*A33);
79a2fa7910SValeria Barra       qdata[i+Q*3] = qw * (A21*A31 + A22*A32 + A23*A33);
80a2fa7910SValeria Barra       qdata[i+Q*4] = qw * (A11*A31 + A12*A32 + A13*A33);
81a2fa7910SValeria Barra       qdata[i+Q*5] = qw * (A11*A21 + A12*A22 + A13*A23);
824d537eeaSYohann     }
834d537eeaSYohann     break;
844d537eeaSYohann   }
854d537eeaSYohann   return 0;
864d537eeaSYohann }
874d537eeaSYohann 
884d537eeaSYohann /// libCEED Q-function for applying a diff operator
894d537eeaSYohann CEED_QFUNCTION(f_apply_diff)(void *ctx, const CeedInt Q,
904d537eeaSYohann                              const CeedScalar *const *in, CeedScalar *const *out) {
914d537eeaSYohann   BuildContext *bc = (BuildContext *)ctx;
924d537eeaSYohann   // in[0], out[0] have shape [dim, nc=1, Q]
93a2fa7910SValeria Barra   const CeedScalar *ug = in[0], *qdata = in[1];
944d537eeaSYohann   CeedScalar *vg = out[0];
95ee07ded2SValeria Barra 
964d537eeaSYohann   switch (bc->dim) {
974d537eeaSYohann   case 1:
98ee07ded2SValeria Barra     // Quadrature Point Loop
99ee07ded2SValeria Barra     CeedPragmaSIMD
1004d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
101a2fa7910SValeria Barra       vg[i] = ug[i] * qdata[i];
1024d537eeaSYohann     }
1034d537eeaSYohann     break;
1044d537eeaSYohann   case 2:
105ee07ded2SValeria Barra     // Quadrature Point Loop
106ee07ded2SValeria Barra     CeedPragmaSIMD
1074d537eeaSYohann     for (CeedInt i=0; i<Q; i++) {
1084d537eeaSYohann       const CeedScalar ug0 = ug[i+Q*0];
1094d537eeaSYohann       const CeedScalar ug1 = ug[i+Q*1];
110a2fa7910SValeria Barra       vg[i+Q*0] = qdata[i+Q*0]*ug0 + qdata[i+Q*2]*ug1;
111a2fa7910SValeria Barra       vg[i+Q*1] = qdata[i+Q*2]*ug0 + qdata[i+Q*1]*ug1;
1124d537eeaSYohann     }
1134d537eeaSYohann     break;
1144d537eeaSYohann   case 3:
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];
1204d537eeaSYohann       const CeedScalar ug2 = ug[i+Q*2];
121a2fa7910SValeria Barra       vg[i+Q*0] = qdata[i+Q*0]*ug0 + qdata[i+Q*5]*ug1 + qdata[i+Q*4]*ug2;
122a2fa7910SValeria Barra       vg[i+Q*1] = qdata[i+Q*5]*ug0 + qdata[i+Q*1]*ug1 + qdata[i+Q*3]*ug2;
123a2fa7910SValeria Barra       vg[i+Q*2] = qdata[i+Q*4]*ug0 + qdata[i+Q*3]*ug1 + qdata[i+Q*2]*ug2;
1244d537eeaSYohann     }
1254d537eeaSYohann     break;
1264d537eeaSYohann   }
1274d537eeaSYohann   return 0;
1284d537eeaSYohann }
129288c0443SJeremy L Thompson 
130a57ca787Svaleriabarra #endif // bp3_h
131