xref: /petsc/src/dm/dt/fe/impls/opencl/feopencl.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
37be5e748SToby Isaac #if defined(PETSC_HAVE_OPENCL)
420cf1dd8SToby Isaac 
52b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_OpenCL(PetscFE fem)
620cf1dd8SToby Isaac {
720cf1dd8SToby Isaac   PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
820cf1dd8SToby Isaac 
920cf1dd8SToby Isaac   PetscFunctionBegin;
105f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseCommandQueue(ocl->queue_id));
1120cf1dd8SToby Isaac   ocl->queue_id = 0;
125f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseContext(ocl->ctx_id));
1320cf1dd8SToby Isaac   ocl->ctx_id = 0;
145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ocl));
1520cf1dd8SToby Isaac   PetscFunctionReturn(0);
1620cf1dd8SToby Isaac }
1720cf1dd8SToby Isaac 
185f80ce2aSJacob Faibussowitsch #define CHKERRSTR(err) do {CHKERRQ(err); string_tail += count; PetscCheck(string_tail != end_of_buffer,PETSC_COMM_SELF, PETSC_ERR_PLIB,"Buffer overflow");} while (0)
1920cf1dd8SToby Isaac enum {LAPLACIAN = 0, ELASTICITY = 1};
2020cf1dd8SToby Isaac 
2120cf1dd8SToby Isaac /* NOTE: This is now broken for vector problems. Must redo loops to respect vector basis elements */
2220cf1dd8SToby Isaac /* dim     Number of spatial dimensions:          2                   */
2320cf1dd8SToby Isaac /* N_b     Number of basis functions:             generated           */
2420cf1dd8SToby Isaac /* N_{bt}  Number of total basis functions:       N_b * N_{comp}      */
2520cf1dd8SToby Isaac /* N_q     Number of quadrature points:           generated           */
2620cf1dd8SToby Isaac /* N_{bs}  Number of block cells                  LCM(N_b, N_q)       */
2720cf1dd8SToby Isaac /* N_{bst} Number of block cell components        LCM(N_{bt}, N_q)    */
2820cf1dd8SToby Isaac /* N_{bl}  Number of concurrent blocks            generated           */
2920cf1dd8SToby Isaac /* N_t     Number of threads:                     N_{bl} * N_{bs}     */
3020cf1dd8SToby Isaac /* N_{cbc} Number of concurrent basis      cells: N_{bl} * N_q        */
3120cf1dd8SToby Isaac /* N_{cqc} Number of concurrent quadrature cells: N_{bl} * N_b        */
3220cf1dd8SToby Isaac /* N_{sbc} Number of serial     basis      cells: N_{bs} / N_q        */
3320cf1dd8SToby Isaac /* N_{sqc} Number of serial     quadrature cells: N_{bs} / N_b        */
3420cf1dd8SToby Isaac /* N_{cb}  Number of serial cell batches:         input               */
3520cf1dd8SToby Isaac /* N_c     Number of total cells:                 N_{cb}*N_{t}/N_{comp} */
362b99622eSMatthew G. Knepley static PetscErrorCode PetscFEOpenCLGenerateIntegrationCode(PetscFE fem, char **string_buffer, PetscInt buffer_length, PetscBool useAux, PetscInt N_bl)
3720cf1dd8SToby Isaac {
3820cf1dd8SToby Isaac   PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
3920cf1dd8SToby Isaac   PetscQuadrature q;
4020cf1dd8SToby Isaac   char           *string_tail   = *string_buffer;
4120cf1dd8SToby Isaac   char           *end_of_buffer = *string_buffer + buffer_length;
4220cf1dd8SToby Isaac   char            float_str[]   = "float", double_str[]  = "double";
4320cf1dd8SToby Isaac   char           *numeric_str   = &(float_str[0]);
4420cf1dd8SToby Isaac   PetscInt        op            = ocl->op;
4520cf1dd8SToby Isaac   PetscBool       useField      = PETSC_FALSE;
4620cf1dd8SToby Isaac   PetscBool       useFieldDer   = PETSC_TRUE;
4720cf1dd8SToby Isaac   PetscBool       useFieldAux   = useAux;
4820cf1dd8SToby Isaac   PetscBool       useFieldDerAux= PETSC_FALSE;
4920cf1dd8SToby Isaac   PetscBool       useF0         = PETSC_TRUE;
5020cf1dd8SToby Isaac   PetscBool       useF1         = PETSC_TRUE;
5120cf1dd8SToby Isaac   const PetscReal *points, *weights;
52ef0bb6c7SMatthew G. Knepley   PetscTabulation T;
5320cf1dd8SToby Isaac   PetscInt        dim, qNc, N_b, N_c, N_q, N_t, p, d, b, c;
5420cf1dd8SToby Isaac   size_t          count;
5520cf1dd8SToby Isaac 
5620cf1dd8SToby Isaac   PetscFunctionBegin;
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetSpatialDimension(fem, &dim));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetDimension(fem, &N_b));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetNumComponents(fem, &N_c));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetQuadrature(fem, &q));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(q, NULL, &qNc, &N_q, &points, &weights));
625f80ce2aSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
6320cf1dd8SToby Isaac   N_t  = N_b * N_c * N_q * N_bl;
6420cf1dd8SToby Isaac   /* Enable device extension for double precision */
6520cf1dd8SToby Isaac   if (ocl->realType == PETSC_DOUBLE) {
665f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
6720cf1dd8SToby Isaac "#if defined(cl_khr_fp64)\n"
6820cf1dd8SToby Isaac "#  pragma OPENCL EXTENSION cl_khr_fp64: enable\n"
6920cf1dd8SToby Isaac "#elif defined(cl_amd_fp64)\n"
7020cf1dd8SToby Isaac "#  pragma OPENCL EXTENSION cl_amd_fp64: enable\n"
7120cf1dd8SToby Isaac "#endif\n",
725f80ce2aSJacob Faibussowitsch                                  &count));
7320cf1dd8SToby Isaac     numeric_str  = &(double_str[0]);
7420cf1dd8SToby Isaac   }
7520cf1dd8SToby Isaac   /* Kernel API */
765f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
7720cf1dd8SToby Isaac "\n"
7820cf1dd8SToby Isaac "__kernel void integrateElementQuadrature(int N_cb, __global %s *coefficients, __global %s *coefficientsAux, __global %s *jacobianInverses, __global %s *jacobianDeterminants, __global %s *elemVec)\n"
7920cf1dd8SToby Isaac "{\n",
805f80ce2aSJacob Faibussowitsch                                &count, numeric_str, numeric_str, numeric_str, numeric_str, numeric_str));
8120cf1dd8SToby Isaac   /* Quadrature */
825f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
8320cf1dd8SToby Isaac "  /* Quadrature points\n"
8420cf1dd8SToby Isaac "   - (x1,y1,x2,y2,...) */\n"
8520cf1dd8SToby Isaac "  const %s points[%d] = {\n",
865f80ce2aSJacob Faibussowitsch                                &count, numeric_str, N_q*dim));
8720cf1dd8SToby Isaac   for (p = 0; p < N_q; ++p) {
8820cf1dd8SToby Isaac     for (d = 0; d < dim; ++d) {
895f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, points[p*dim+d]));
9020cf1dd8SToby Isaac     }
9120cf1dd8SToby Isaac   }
925f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count));
935f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
9420cf1dd8SToby Isaac "  /* Quadrature weights\n"
9520cf1dd8SToby Isaac "   - (v1,v2,...) */\n"
9620cf1dd8SToby Isaac "  const %s weights[%d] = {\n",
975f80ce2aSJacob Faibussowitsch                                &count, numeric_str, N_q));
9820cf1dd8SToby Isaac   for (p = 0; p < N_q; ++p) {
995f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, weights[p]));
10020cf1dd8SToby Isaac   }
1015f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count));
10220cf1dd8SToby Isaac   /* Basis Functions */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetCellTabulation(fem, 1, &T));
1045f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
10520cf1dd8SToby Isaac "  /* Nodal basis function evaluations\n"
10620cf1dd8SToby Isaac "    - basis component is fastest varying, the basis function, then point */\n"
10720cf1dd8SToby Isaac "  const %s Basis[%d] = {\n",
1085f80ce2aSJacob Faibussowitsch                                &count, numeric_str, N_q*N_b*N_c));
10920cf1dd8SToby Isaac   for (p = 0; p < N_q; ++p) {
11020cf1dd8SToby Isaac     for (b = 0; b < N_b; ++b) {
11120cf1dd8SToby Isaac       for (c = 0; c < N_c; ++c) {
1125f80ce2aSJacob Faibussowitsch         CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, T->T[0][(p*N_b + b)*N_c + c]));
11320cf1dd8SToby Isaac       }
11420cf1dd8SToby Isaac     }
11520cf1dd8SToby Isaac   }
1165f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count));
1175f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
11820cf1dd8SToby Isaac "\n"
11920cf1dd8SToby Isaac "  /* Nodal basis function derivative evaluations,\n"
12020cf1dd8SToby Isaac "      - derivative direction is fastest varying, then basis component, then basis function, then point */\n"
12120cf1dd8SToby Isaac "  const %s%d BasisDerivatives[%d] = {\n",
1225f80ce2aSJacob Faibussowitsch                             &count, numeric_str, dim, N_q*N_b*N_c));
12320cf1dd8SToby Isaac   for (p = 0; p < N_q; ++p) {
12420cf1dd8SToby Isaac     for (b = 0; b < N_b; ++b) {
12520cf1dd8SToby Isaac       for (c = 0; c < N_c; ++c) {
1265f80ce2aSJacob Faibussowitsch         CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "(%s%d)(", &count, numeric_str, dim));
12720cf1dd8SToby Isaac         for (d = 0; d < dim; ++d) {
12820cf1dd8SToby Isaac           if (d > 0) {
1295f80ce2aSJacob Faibussowitsch             CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, ", %g", &count, T->T[1][((p*N_b + b)*dim + d)*N_c + c]));
13020cf1dd8SToby Isaac           } else {
1315f80ce2aSJacob Faibussowitsch             CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g", &count, T->T[1][((p*N_b + b)*dim + d)*N_c + c]));
13220cf1dd8SToby Isaac           }
13320cf1dd8SToby Isaac         }
1345f80ce2aSJacob Faibussowitsch         CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "),\n", &count));
13520cf1dd8SToby Isaac       }
13620cf1dd8SToby Isaac     }
13720cf1dd8SToby Isaac   }
1385f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count));
13920cf1dd8SToby Isaac   /* Sizes */
1405f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
14120cf1dd8SToby Isaac "  const int dim    = %d;                           // The spatial dimension\n"
14220cf1dd8SToby Isaac "  const int N_bl   = %d;                           // The number of concurrent blocks\n"
14320cf1dd8SToby Isaac "  const int N_b    = %d;                           // The number of basis functions\n"
14420cf1dd8SToby Isaac "  const int N_comp = %d;                           // The number of basis function components\n"
14520cf1dd8SToby Isaac "  const int N_bt   = N_b*N_comp;                    // The total number of scalar basis functions\n"
14620cf1dd8SToby Isaac "  const int N_q    = %d;                           // The number of quadrature points\n"
14720cf1dd8SToby Isaac "  const int N_bst  = N_bt*N_q;                      // The block size, LCM(N_b*N_comp, N_q), Notice that a block is not processed simultaneously\n"
14820cf1dd8SToby Isaac "  const int N_t    = N_bst*N_bl;                    // The number of threads, N_bst * N_bl\n"
14920cf1dd8SToby Isaac "  const int N_bc   = N_t/N_comp;                    // The number of cells per batch (N_b*N_q*N_bl)\n"
15020cf1dd8SToby Isaac "  const int N_sbc  = N_bst / (N_q * N_comp);\n"
15120cf1dd8SToby Isaac "  const int N_sqc  = N_bst / N_bt;\n"
15220cf1dd8SToby Isaac "  /*const int N_c    = N_cb * N_bc;*/\n"
15320cf1dd8SToby Isaac "\n"
15420cf1dd8SToby Isaac "  /* Calculated indices */\n"
15520cf1dd8SToby Isaac "  /*const int tidx    = get_local_id(0) + get_local_size(0)*get_local_id(1);*/\n"
15620cf1dd8SToby Isaac "  const int tidx    = get_local_id(0);\n"
15720cf1dd8SToby Isaac "  const int blidx   = tidx / N_bst;                  // Block number for this thread\n"
15820cf1dd8SToby Isaac "  const int bidx    = tidx %% N_bt;                   // Basis function mapped to this thread\n"
15920cf1dd8SToby Isaac "  const int cidx    = tidx %% N_comp;                 // Basis component mapped to this thread\n"
16020cf1dd8SToby Isaac "  const int qidx    = tidx %% N_q;                    // Quadrature point mapped to this thread\n"
16120cf1dd8SToby Isaac "  const int blbidx  = tidx %% N_q + blidx*N_q;        // Cell mapped to this thread in the basis phase\n"
16220cf1dd8SToby Isaac "  const int blqidx  = tidx %% N_b + blidx*N_b;        // Cell mapped to this thread in the quadrature phase\n"
16320cf1dd8SToby Isaac "  const int gidx    = get_group_id(1)*get_num_groups(0) + get_group_id(0);\n"
16420cf1dd8SToby Isaac "  const int Goffset = gidx*N_cb*N_bc;\n",
1655f80ce2aSJacob Faibussowitsch                                &count, dim, N_bl, N_b, N_c, N_q));
16620cf1dd8SToby Isaac   /* Local memory */
1675f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
16820cf1dd8SToby Isaac "\n"
16920cf1dd8SToby Isaac "  /* Quadrature data */\n"
17020cf1dd8SToby Isaac "  %s                w;                   // $w_q$, Quadrature weight at $x_q$\n"
17120cf1dd8SToby Isaac "  __local %s         phi_i[%d];    //[N_bt*N_q];  // $\\phi_i(x_q)$, Value of the basis function $i$ at $x_q$\n"
17220cf1dd8SToby Isaac "  __local %s%d       phiDer_i[%d]; //[N_bt*N_q];  // $\\frac{\\partial\\phi_i(x_q)}{\\partial x_d}$, Value of the derivative of basis function $i$ in direction $x_d$ at $x_q$\n"
17320cf1dd8SToby Isaac "  /* Geometric data */\n"
17420cf1dd8SToby Isaac "  __local %s        detJ[%d]; //[N_t];           // $|J(x_q)|$, Jacobian determinant at $x_q$\n"
17520cf1dd8SToby Isaac "  __local %s        invJ[%d];//[N_t*dim*dim];   // $J^{-1}(x_q)$, Jacobian inverse at $x_q$\n",
17620cf1dd8SToby Isaac                                &count, numeric_str, numeric_str, N_b*N_c*N_q, numeric_str, dim, N_b*N_c*N_q, numeric_str, N_t,
1775f80ce2aSJacob Faibussowitsch                                numeric_str, N_t*dim*dim, numeric_str, N_t*N_b*N_c));
1785f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
17920cf1dd8SToby Isaac "  /* FEM data */\n"
18020cf1dd8SToby Isaac "  __local %s        u_i[%d]; //[N_t*N_bt];       // Coefficients $u_i$ of the field $u|_{\\mathcal{T}} = \\sum_i u_i \\phi_i$\n",
1815f80ce2aSJacob Faibussowitsch                                &count, numeric_str, N_t*N_b*N_c));
18220cf1dd8SToby Isaac   if (useAux) {
1835f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
18420cf1dd8SToby Isaac "  __local %s        a_i[%d]; //[N_t];            // Coefficients $a_i$ of the auxiliary field $a|_{\\mathcal{T}} = \\sum_i a_i \\phi^R_i$\n",
1855f80ce2aSJacob Faibussowitsch                                  &count, numeric_str, N_t));
18620cf1dd8SToby Isaac   }
18720cf1dd8SToby Isaac   if (useF0) {
1885f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
18920cf1dd8SToby Isaac "  /* Intermediate calculations */\n"
19020cf1dd8SToby Isaac "  __local %s         f_0[%d]; //[N_t*N_sqc];      // $f_0(u(x_q), \\nabla u(x_q)) |J(x_q)| w_q$\n",
1915f80ce2aSJacob Faibussowitsch                                  &count, numeric_str, N_t*N_q));
19220cf1dd8SToby Isaac   }
19320cf1dd8SToby Isaac   if (useF1) {
1945f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
19520cf1dd8SToby Isaac "  __local %s%d       f_1[%d]; //[N_t*N_sqc];      // $f_1(u(x_q), \\nabla u(x_q)) |J(x_q)| w_q$\n",
1965f80ce2aSJacob Faibussowitsch                                  &count, numeric_str, dim, N_t*N_q));
19720cf1dd8SToby Isaac   }
19820cf1dd8SToby Isaac   /* TODO: If using elasticity, put in mu/lambda coefficients */
1995f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
20020cf1dd8SToby Isaac "  /* Output data */\n"
20120cf1dd8SToby Isaac "  %s                e_i;                 // Coefficient $e_i$ of the residual\n\n",
2025f80ce2aSJacob Faibussowitsch                                &count, numeric_str));
20320cf1dd8SToby Isaac   /* One-time loads */
2045f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
20520cf1dd8SToby Isaac "  /* These should be generated inline */\n"
20620cf1dd8SToby Isaac "  /* Load quadrature weights */\n"
20720cf1dd8SToby Isaac "  w = weights[qidx];\n"
20820cf1dd8SToby Isaac "  /* Load basis tabulation \\phi_i for this cell */\n"
20920cf1dd8SToby Isaac "  if (tidx < N_bt*N_q) {\n"
21020cf1dd8SToby Isaac "    phi_i[tidx]    = Basis[tidx];\n"
21120cf1dd8SToby Isaac "    phiDer_i[tidx] = BasisDerivatives[tidx];\n"
21220cf1dd8SToby Isaac "  }\n\n",
2135f80ce2aSJacob Faibussowitsch                                &count));
21420cf1dd8SToby Isaac   /* Batch loads */
2155f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
21620cf1dd8SToby Isaac "  for (int batch = 0; batch < N_cb; ++batch) {\n"
21720cf1dd8SToby Isaac "    /* Load geometry */\n"
21820cf1dd8SToby Isaac "    detJ[tidx] = jacobianDeterminants[Goffset+batch*N_bc+tidx];\n"
21920cf1dd8SToby Isaac "    for (int n = 0; n < dim*dim; ++n) {\n"
22020cf1dd8SToby Isaac "      const int offset = n*N_t;\n"
22120cf1dd8SToby Isaac "      invJ[offset+tidx] = jacobianInverses[(Goffset+batch*N_bc)*dim*dim+offset+tidx];\n"
22220cf1dd8SToby Isaac "    }\n"
22320cf1dd8SToby Isaac "    /* Load coefficients u_i for this cell */\n"
22420cf1dd8SToby Isaac "    for (int n = 0; n < N_bt; ++n) {\n"
22520cf1dd8SToby Isaac "      const int offset = n*N_t;\n"
22620cf1dd8SToby Isaac "      u_i[offset+tidx] = coefficients[(Goffset*N_bt)+batch*N_t*N_b+offset+tidx];\n"
22720cf1dd8SToby Isaac "    }\n",
2285f80ce2aSJacob Faibussowitsch                                &count));
22920cf1dd8SToby Isaac   if (useAux) {
2305f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
23120cf1dd8SToby Isaac "    /* Load coefficients a_i for this cell */\n"
23220cf1dd8SToby Isaac "    /* TODO: This should not be N_t here, it should be N_bc*N_comp_aux */\n"
23320cf1dd8SToby Isaac "    a_i[tidx] = coefficientsAux[Goffset+batch*N_t+tidx];\n",
2345f80ce2aSJacob Faibussowitsch                                  &count));
23520cf1dd8SToby Isaac   }
23620cf1dd8SToby Isaac   /* Quadrature phase */
2375f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
23820cf1dd8SToby Isaac "    barrier(CLK_LOCAL_MEM_FENCE);\n"
23920cf1dd8SToby Isaac "\n"
24020cf1dd8SToby Isaac "    /* Map coefficients to values at quadrature points */\n"
24120cf1dd8SToby Isaac "    for (int c = 0; c < N_sqc; ++c) {\n"
24220cf1dd8SToby Isaac "      const int cell          = c*N_bl*N_b + blqidx;\n"
24320cf1dd8SToby Isaac "      const int fidx          = (cell*N_q + qidx)*N_comp + cidx;\n",
2445f80ce2aSJacob Faibussowitsch                                &count));
24520cf1dd8SToby Isaac   if (useField) {
2465f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
24720cf1dd8SToby Isaac "      %s  u[%d]; //[N_comp];     // $u(x_q)$, Value of the field at $x_q$\n",
2485f80ce2aSJacob Faibussowitsch                                  &count, numeric_str, N_c));
24920cf1dd8SToby Isaac   }
25020cf1dd8SToby Isaac   if (useFieldDer) {
2515f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
25220cf1dd8SToby Isaac "      %s%d   gradU[%d]; //[N_comp]; // $\\nabla u(x_q)$, Value of the field gradient at $x_q$\n",
2535f80ce2aSJacob Faibussowitsch                                  &count, numeric_str, dim, N_c));
25420cf1dd8SToby Isaac   }
25520cf1dd8SToby Isaac   if (useFieldAux) {
2565f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
25720cf1dd8SToby Isaac "      %s  a[%d]; //[1];     // $a(x_q)$, Value of the auxiliary fields at $x_q$\n",
2585f80ce2aSJacob Faibussowitsch                                  &count, numeric_str, 1));
25920cf1dd8SToby Isaac   }
26020cf1dd8SToby Isaac   if (useFieldDerAux) {
2615f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
26220cf1dd8SToby Isaac "      %s%d   gradA[%d]; //[1]; // $\\nabla a(x_q)$, Value of the auxiliary field gradient at $x_q$\n",
2635f80ce2aSJacob Faibussowitsch                                  &count, numeric_str, dim, 1));
26420cf1dd8SToby Isaac   }
2655f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
26620cf1dd8SToby Isaac "\n"
26720cf1dd8SToby Isaac "      for (int comp = 0; comp < N_comp; ++comp) {\n",
2685f80ce2aSJacob Faibussowitsch                                &count));
2695f80ce2aSJacob Faibussowitsch   if (useField) CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "        u[comp] = 0.0;\n", &count));
27020cf1dd8SToby Isaac   if (useFieldDer) {
27120cf1dd8SToby Isaac     switch (dim) {
27220cf1dd8SToby Isaac     case 1:
2735f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "        gradU[comp].x = 0.0;\n", &count));break;
27420cf1dd8SToby Isaac     case 2:
2755f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "        gradU[comp].x = 0.0; gradU[comp].y = 0.0;\n", &count));break;
27620cf1dd8SToby Isaac     case 3:
2775f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "        gradU[comp].x = 0.0; gradU[comp].y = 0.0; gradU[comp].z = 0.0;\n", &count));break;
27820cf1dd8SToby Isaac     }
27920cf1dd8SToby Isaac   }
2805f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
28120cf1dd8SToby Isaac "      }\n",
2825f80ce2aSJacob Faibussowitsch                                &count));
28320cf1dd8SToby Isaac   if (useFieldAux) {
2845f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      a[0] = 0.0;\n", &count));
28520cf1dd8SToby Isaac   }
28620cf1dd8SToby Isaac   if (useFieldDerAux) {
28720cf1dd8SToby Isaac     switch (dim) {
28820cf1dd8SToby Isaac     case 1:
2895f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      gradA[0].x = 0.0;\n", &count));break;
29020cf1dd8SToby Isaac     case 2:
2915f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      gradA[0].x = 0.0; gradA[0].y = 0.0;\n", &count));break;
29220cf1dd8SToby Isaac     case 3:
2935f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      gradA[0].x = 0.0; gradA[0].y = 0.0; gradA[0].z = 0.0;\n", &count));break;
29420cf1dd8SToby Isaac     }
29520cf1dd8SToby Isaac   }
2965f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
29720cf1dd8SToby Isaac "      /* Get field and derivatives at this quadrature point */\n"
29820cf1dd8SToby Isaac "      for (int i = 0; i < N_b; ++i) {\n"
29920cf1dd8SToby Isaac "        for (int comp = 0; comp < N_comp; ++comp) {\n"
30020cf1dd8SToby Isaac "          const int b    = i*N_comp+comp;\n"
30120cf1dd8SToby Isaac "          const int pidx = qidx*N_bt + b;\n"
30220cf1dd8SToby Isaac "          const int uidx = cell*N_bt + b;\n"
30320cf1dd8SToby Isaac "          %s%d   realSpaceDer;\n\n",
3045f80ce2aSJacob Faibussowitsch                                &count, numeric_str, dim));
3055f80ce2aSJacob Faibussowitsch   if (useField) CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,"          u[comp] += u_i[uidx]*phi_i[pidx];\n", &count));
30620cf1dd8SToby Isaac   if (useFieldDer) {
30720cf1dd8SToby Isaac     switch (dim) {
30820cf1dd8SToby Isaac     case 2:
3095f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
31020cf1dd8SToby Isaac "          realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y;\n"
31120cf1dd8SToby Isaac "          gradU[comp].x += u_i[uidx]*realSpaceDer.x;\n"
31220cf1dd8SToby Isaac "          realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y;\n"
31320cf1dd8SToby Isaac "          gradU[comp].y += u_i[uidx]*realSpaceDer.y;\n",
3145f80ce2aSJacob Faibussowitsch                                    &count));break;
31520cf1dd8SToby Isaac     case 3:
3165f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
31720cf1dd8SToby Isaac "          realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+0]*phiDer_i[pidx].z;\n"
31820cf1dd8SToby Isaac "          gradU[comp].x += u_i[uidx]*realSpaceDer.x;\n"
31920cf1dd8SToby Isaac "          realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+1]*phiDer_i[pidx].z;\n"
32020cf1dd8SToby Isaac "          gradU[comp].y += u_i[uidx]*realSpaceDer.y;\n"
32120cf1dd8SToby Isaac "          realSpaceDer.z = invJ[cell*dim*dim+0*dim+2]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+2]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+2]*phiDer_i[pidx].z;\n"
32220cf1dd8SToby Isaac "          gradU[comp].z += u_i[uidx]*realSpaceDer.z;\n",
3235f80ce2aSJacob Faibussowitsch                                    &count));break;
32420cf1dd8SToby Isaac     }
32520cf1dd8SToby Isaac   }
3265f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
32720cf1dd8SToby Isaac "        }\n"
32820cf1dd8SToby Isaac "      }\n",
3295f80ce2aSJacob Faibussowitsch                                &count));
33020cf1dd8SToby Isaac   if (useFieldAux) {
3315f80ce2aSJacob Faibussowitsch     CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,"          a[0] += a_i[cell];\n", &count));
33220cf1dd8SToby Isaac   }
33320cf1dd8SToby Isaac   /* Calculate residual at quadrature points: Should be generated by an weak form egine */
3345f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
33520cf1dd8SToby Isaac "      /* Process values at quadrature points */\n",
3365f80ce2aSJacob Faibussowitsch                                &count));
33720cf1dd8SToby Isaac   switch (op) {
33820cf1dd8SToby Isaac   case LAPLACIAN:
3395f80ce2aSJacob Faibussowitsch     if (useF0) {CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      f_0[fidx] = 4.0;\n", &count));}
34020cf1dd8SToby Isaac     if (useF1) {
3415f80ce2aSJacob Faibussowitsch       if (useAux) CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      f_1[fidx] = a[0]*gradU[cidx];\n", &count));
3425f80ce2aSJacob Faibussowitsch       else        CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      f_1[fidx] = gradU[cidx];\n", &count));
34320cf1dd8SToby Isaac     }
34420cf1dd8SToby Isaac     break;
34520cf1dd8SToby Isaac   case ELASTICITY:
3465f80ce2aSJacob Faibussowitsch     if (useF0) CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "      f_0[fidx] = 4.0;\n", &count));
34720cf1dd8SToby Isaac     if (useF1) {
34820cf1dd8SToby Isaac     switch (dim) {
34920cf1dd8SToby Isaac     case 2:
3505f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
35120cf1dd8SToby Isaac "      switch (cidx) {\n"
35220cf1dd8SToby Isaac "      case 0:\n"
35320cf1dd8SToby Isaac "        f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[0].x + gradU[0].x);\n"
35420cf1dd8SToby Isaac "        f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[0].y + gradU[1].x);\n"
35520cf1dd8SToby Isaac "        break;\n"
35620cf1dd8SToby Isaac "      case 1:\n"
35720cf1dd8SToby Isaac "        f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[1].x + gradU[0].y);\n"
35820cf1dd8SToby Isaac "        f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[1].y + gradU[1].y);\n"
35920cf1dd8SToby Isaac "      }\n",
3605f80ce2aSJacob Faibussowitsch                                    &count));break;
36120cf1dd8SToby Isaac     case 3:
3625f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
36320cf1dd8SToby Isaac "      switch (cidx) {\n"
36420cf1dd8SToby Isaac "      case 0:\n"
36520cf1dd8SToby Isaac "        f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].x + gradU[0].x);\n"
36620cf1dd8SToby Isaac "        f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].y + gradU[1].x);\n"
36720cf1dd8SToby Isaac "        f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].z + gradU[2].x);\n"
36820cf1dd8SToby Isaac "        break;\n"
36920cf1dd8SToby Isaac "      case 1:\n"
37020cf1dd8SToby Isaac "        f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].x + gradU[0].y);\n"
37120cf1dd8SToby Isaac "        f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].y + gradU[1].y);\n"
37220cf1dd8SToby Isaac "        f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].y + gradU[2].y);\n"
37320cf1dd8SToby Isaac "        break;\n"
37420cf1dd8SToby Isaac "      case 2:\n"
37520cf1dd8SToby Isaac "        f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].x + gradU[0].z);\n"
37620cf1dd8SToby Isaac "        f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].y + gradU[1].z);\n"
37720cf1dd8SToby Isaac "        f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].y + gradU[2].z);\n"
37820cf1dd8SToby Isaac "      }\n",
3795f80ce2aSJacob Faibussowitsch                                    &count));break;
38020cf1dd8SToby Isaac     }}
38120cf1dd8SToby Isaac     break;
38220cf1dd8SToby Isaac   default:
38398921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "PDE operator %d is not supported", op);
38420cf1dd8SToby Isaac   }
3855f80ce2aSJacob Faibussowitsch   if (useF0) CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,"      f_0[fidx] *= detJ[cell]*w;\n", &count));
38620cf1dd8SToby Isaac   if (useF1) {
38720cf1dd8SToby Isaac     switch (dim) {
38820cf1dd8SToby Isaac     case 1:
3895f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,"      f_1[fidx].x *= detJ[cell]*w;\n", &count));break;
39020cf1dd8SToby Isaac     case 2:
3915f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,"      f_1[fidx].x *= detJ[cell]*w; f_1[fidx].y *= detJ[cell]*w;\n", &count));break;
39220cf1dd8SToby Isaac     case 3:
3935f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,"      f_1[fidx].x *= detJ[cell]*w; f_1[fidx].y *= detJ[cell]*w; f_1[fidx].z *= detJ[cell]*w;\n", &count));break;
39420cf1dd8SToby Isaac     }
39520cf1dd8SToby Isaac   }
39620cf1dd8SToby Isaac   /* Thread transpose */
3975f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
39820cf1dd8SToby Isaac "    }\n\n"
39920cf1dd8SToby Isaac "    /* ==== TRANSPOSE THREADS ==== */\n"
40020cf1dd8SToby Isaac "    barrier(CLK_LOCAL_MEM_FENCE);\n\n",
4015f80ce2aSJacob Faibussowitsch                                &count));
40220cf1dd8SToby Isaac   /* Basis phase */
4035f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
40420cf1dd8SToby Isaac "    /* Map values at quadrature points to coefficients */\n"
40520cf1dd8SToby Isaac "    for (int c = 0; c < N_sbc; ++c) {\n"
40620cf1dd8SToby Isaac "      const int cell = c*N_bl*N_q + blbidx; /* Cell number in batch */\n"
40720cf1dd8SToby Isaac "\n"
40820cf1dd8SToby Isaac "      e_i = 0.0;\n"
40920cf1dd8SToby Isaac "      for (int q = 0; q < N_q; ++q) {\n"
41020cf1dd8SToby Isaac "        const int pidx = q*N_bt + bidx;\n"
41120cf1dd8SToby Isaac "        const int fidx = (cell*N_q + q)*N_comp + cidx;\n"
41220cf1dd8SToby Isaac "        %s%d   realSpaceDer;\n\n",
4135f80ce2aSJacob Faibussowitsch                                &count, numeric_str, dim));
41420cf1dd8SToby Isaac 
4155f80ce2aSJacob Faibussowitsch   if (useF0) CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,"        e_i += phi_i[pidx]*f_0[fidx];\n", &count));
41620cf1dd8SToby Isaac   if (useF1) {
41720cf1dd8SToby Isaac     switch (dim) {
41820cf1dd8SToby Isaac     case 2:
4195f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
42020cf1dd8SToby Isaac "        realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y;\n"
42120cf1dd8SToby Isaac "        e_i           += realSpaceDer.x*f_1[fidx].x;\n"
42220cf1dd8SToby Isaac "        realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y;\n"
42320cf1dd8SToby Isaac "        e_i           += realSpaceDer.y*f_1[fidx].y;\n",
4245f80ce2aSJacob Faibussowitsch                                    &count));break;
42520cf1dd8SToby Isaac     case 3:
4265f80ce2aSJacob Faibussowitsch       CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
42720cf1dd8SToby Isaac "        realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+0]*phiDer_i[pidx].z;\n"
42820cf1dd8SToby Isaac "        e_i           += realSpaceDer.x*f_1[fidx].x;\n"
42920cf1dd8SToby Isaac "        realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+1]*phiDer_i[pidx].z;\n"
43020cf1dd8SToby Isaac "        e_i           += realSpaceDer.y*f_1[fidx].y;\n"
43120cf1dd8SToby Isaac "        realSpaceDer.z = invJ[cell*dim*dim+0*dim+2]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+2]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+2]*phiDer_i[pidx].z;\n"
43220cf1dd8SToby Isaac "        e_i           += realSpaceDer.z*f_1[fidx].z;\n",
4335f80ce2aSJacob Faibussowitsch                                    &count));break;
43420cf1dd8SToby Isaac     }
43520cf1dd8SToby Isaac   }
4365f80ce2aSJacob Faibussowitsch   CHKERRSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
43720cf1dd8SToby Isaac "      }\n"
43820cf1dd8SToby Isaac "      /* Write element vector for N_{cbc} cells at a time */\n"
43920cf1dd8SToby Isaac "      elemVec[(Goffset + batch*N_bc + c*N_bl*N_q)*N_bt + tidx] = e_i;\n"
44020cf1dd8SToby Isaac "    }\n"
44120cf1dd8SToby Isaac "    /* ==== Could do one write per batch ==== */\n"
44220cf1dd8SToby Isaac "  }\n"
44320cf1dd8SToby Isaac "  return;\n"
44420cf1dd8SToby Isaac "}\n",
4455f80ce2aSJacob Faibussowitsch                                &count));
44620cf1dd8SToby Isaac   PetscFunctionReturn(0);
44720cf1dd8SToby Isaac }
44820cf1dd8SToby Isaac 
4492b99622eSMatthew G. Knepley static PetscErrorCode PetscFEOpenCLGetIntegrationKernel(PetscFE fem, PetscBool useAux, cl_program *ocl_prog, cl_kernel *ocl_kernel)
45020cf1dd8SToby Isaac {
45120cf1dd8SToby Isaac   PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
45220cf1dd8SToby Isaac   PetscInt        dim, N_bl;
45320cf1dd8SToby Isaac   PetscBool       flg;
45420cf1dd8SToby Isaac   char           *buffer;
45520cf1dd8SToby Isaac   size_t          len;
45620cf1dd8SToby Isaac   char            errMsg[8192];
4572da392ccSBarry Smith   cl_int          err;
45820cf1dd8SToby Isaac   PetscErrorCode  ierr;
45920cf1dd8SToby Isaac 
46020cf1dd8SToby Isaac   PetscFunctionBegin;
4615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetSpatialDimension(fem, &dim));
4625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(8192, &buffer));
4635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetTileSizes(fem, NULL, &N_bl, NULL, NULL));
4645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEOpenCLGenerateIntegrationCode(fem, &buffer, 8192, useAux, N_bl));
4655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(((PetscObject)fem)->options,((PetscObject)fem)->prefix, "-petscfe_opencl_kernel_print", &flg));
4665f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) fem), "OpenCL FE Integration Kernel:\n%s\n", buffer));
4675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlen(buffer,&len));
4682da392ccSBarry Smith   *ocl_prog = clCreateProgramWithSource(ocl->ctx_id, 1, (const char **) &buffer, &len, &err);CHKERRQ(err);
4692da392ccSBarry Smith   err = clBuildProgram(*ocl_prog, 0, NULL, NULL, NULL, NULL);
4702da392ccSBarry Smith   if (err != CL_SUCCESS) {
4712f613bf5SBarry Smith     err = clGetProgramBuildInfo(*ocl_prog, ocl->dev_id, CL_PROGRAM_BUILD_LOG, 8192*sizeof(char), &errMsg, NULL);
47298921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Build failed! Log:\n %s", errMsg);
47320cf1dd8SToby Isaac   }
4745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(buffer));
4752f613bf5SBarry Smith   *ocl_kernel = clCreateKernel(*ocl_prog, "integrateElementQuadrature", &ierr);
47620cf1dd8SToby Isaac   PetscFunctionReturn(0);
47720cf1dd8SToby Isaac }
47820cf1dd8SToby Isaac 
4792b99622eSMatthew G. Knepley static PetscErrorCode PetscFEOpenCLCalculateGrid(PetscFE fem, PetscInt N, PetscInt blockSize, size_t *x, size_t *y, size_t *z)
48020cf1dd8SToby Isaac {
48120cf1dd8SToby Isaac   const PetscInt Nblocks = N/blockSize;
48220cf1dd8SToby Isaac 
48320cf1dd8SToby Isaac   PetscFunctionBegin;
4845f80ce2aSJacob Faibussowitsch   PetscCheck(!(N % blockSize),PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid block size %d for %d elements", blockSize, N);
48520cf1dd8SToby Isaac   *z = 1;
486623c9f23SSatish Balay   *y = 1;
48720cf1dd8SToby Isaac   for (*x = (size_t) (PetscSqrtReal(Nblocks) + 0.5); *x > 0; --*x) {
48820cf1dd8SToby Isaac     *y = Nblocks / *x;
489526996e7SStefano Zampini     if (*x * *y == (size_t)Nblocks) break;
49020cf1dd8SToby Isaac   }
4915f80ce2aSJacob Faibussowitsch   PetscCheck(*x * *y == (size_t)Nblocks,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Could not find partition for %D with block size %D", N, blockSize);
49220cf1dd8SToby Isaac   PetscFunctionReturn(0);
49320cf1dd8SToby Isaac }
49420cf1dd8SToby Isaac 
4952b99622eSMatthew G. Knepley static PetscErrorCode PetscFEOpenCLLogResidual(PetscFE fem, PetscLogDouble time, PetscLogDouble flops)
49620cf1dd8SToby Isaac {
49720cf1dd8SToby Isaac   PetscFE_OpenCL   *ocl = (PetscFE_OpenCL *) fem->data;
49820cf1dd8SToby Isaac   PetscStageLog     stageLog;
49920cf1dd8SToby Isaac   PetscEventPerfLog eventLog = NULL;
500afbfd3a7SStefano Zampini   int               stage;
50120cf1dd8SToby Isaac 
50220cf1dd8SToby Isaac   PetscFunctionBegin;
5035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogGetStageLog(&stageLog));
5045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStageLogGetCurrent(stageLog, &stage));
5055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStageLogGetEventPerfLog(stageLog, stage, &eventLog));
50620cf1dd8SToby Isaac     /* Log performance info */
50720cf1dd8SToby Isaac   eventLog->eventInfo[ocl->residualEvent].count++;
50820cf1dd8SToby Isaac   eventLog->eventInfo[ocl->residualEvent].time  += time;
50920cf1dd8SToby Isaac   eventLog->eventInfo[ocl->residualEvent].flops += flops;
51020cf1dd8SToby Isaac   PetscFunctionReturn(0);
51120cf1dd8SToby Isaac }
51220cf1dd8SToby Isaac 
51306ad1575SMatthew G. Knepley static PetscErrorCode PetscFEIntegrateResidual_OpenCL(PetscDS prob, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
51420cf1dd8SToby Isaac                                                       const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
51520cf1dd8SToby Isaac {
51620cf1dd8SToby Isaac   /* Nbc = batchSize */
517849189efSMatthew G. Knepley   PetscFE           fem;
518849189efSMatthew G. Knepley   PetscFE_OpenCL   *ocl;
51920cf1dd8SToby Isaac   PetscPointFunc    f0_func;
52020cf1dd8SToby Isaac   PetscPointFunc    f1_func;
52120cf1dd8SToby Isaac   PetscQuadrature   q;
52220cf1dd8SToby Isaac   PetscInt          dim, qNc;
52320cf1dd8SToby Isaac   PetscInt          N_b;    /* The number of basis functions */
52420cf1dd8SToby Isaac   PetscInt          N_comp; /* The number of basis function components */
52520cf1dd8SToby Isaac   PetscInt          N_bt;   /* The total number of scalar basis functions */
52620cf1dd8SToby Isaac   PetscInt          N_q;    /* The number of quadrature points */
52720cf1dd8SToby Isaac   PetscInt          N_bst;  /* The block size, LCM(N_bt, N_q), Notice that a block is not process simultaneously */
52820cf1dd8SToby Isaac   PetscInt          N_t;    /* The number of threads, N_bst * N_bl */
52920cf1dd8SToby Isaac   PetscInt          N_bl;   /* The number of blocks */
53020cf1dd8SToby Isaac   PetscInt          N_bc;   /* The batch size, N_bl*N_q*N_b */
53120cf1dd8SToby Isaac   PetscInt          N_cb;   /* The number of batches */
5326528b96dSMatthew G. Knepley   const PetscInt    field = key.field;
53320cf1dd8SToby Isaac   PetscInt          numFlops, f0Flops = 0, f1Flops = 0;
53420cf1dd8SToby Isaac   PetscBool         useAux      = probAux ? PETSC_TRUE : PETSC_FALSE;
53520cf1dd8SToby Isaac   PetscBool         useField    = PETSC_FALSE;
53620cf1dd8SToby Isaac   PetscBool         useFieldDer = PETSC_TRUE;
53720cf1dd8SToby Isaac   PetscBool         useF0       = PETSC_TRUE;
53820cf1dd8SToby Isaac   PetscBool         useF1       = PETSC_TRUE;
53920cf1dd8SToby Isaac   /* OpenCL variables */
54020cf1dd8SToby Isaac   cl_program        ocl_prog;
54120cf1dd8SToby Isaac   cl_kernel         ocl_kernel;
54220cf1dd8SToby Isaac   cl_event          ocl_ev;         /* The event for tracking kernel execution */
54320cf1dd8SToby Isaac   cl_ulong          ns_start;       /* Nanoseconds counter on GPU at kernel start */
54420cf1dd8SToby Isaac   cl_ulong          ns_end;         /* Nanoseconds counter on GPU at kernel stop */
54520cf1dd8SToby Isaac   cl_mem            o_jacobianInverses, o_jacobianDeterminants;
54620cf1dd8SToby Isaac   cl_mem            o_coefficients, o_coefficientsAux, o_elemVec;
54720cf1dd8SToby Isaac   float            *f_coeff = NULL, *f_coeffAux = NULL, *f_invJ = NULL, *f_detJ = NULL;
54820cf1dd8SToby Isaac   double           *d_coeff = NULL, *d_coeffAux = NULL, *d_invJ = NULL, *d_detJ = NULL;
54920cf1dd8SToby Isaac   PetscReal        *r_invJ = NULL, *r_detJ = NULL;
55020cf1dd8SToby Isaac   void             *oclCoeff, *oclCoeffAux, *oclInvJ, *oclDetJ;
55120cf1dd8SToby Isaac   size_t            local_work_size[3], global_work_size[3];
55220cf1dd8SToby Isaac   size_t            realSize, x, y, z;
55320cf1dd8SToby Isaac   const PetscReal   *points, *weights;
55420cf1dd8SToby Isaac   PetscErrorCode    ierr;
55520cf1dd8SToby Isaac 
55620cf1dd8SToby Isaac   PetscFunctionBegin;
5575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetDiscretization(prob, field, (PetscObject *) &fem));
558849189efSMatthew G. Knepley   ocl  = (PetscFE_OpenCL *) fem->data;
5595f80ce2aSJacob Faibussowitsch   if (!Ne) {CHKERRQ(PetscFEOpenCLLogResidual(fem, 0.0, 0.0)); PetscFunctionReturn(0);}
5605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetSpatialDimension(fem, &dim));
5615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetQuadrature(fem, &q));
5625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(q, NULL, &qNc, &N_q, &points, &weights));
5635f80ce2aSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components", qNc);
5645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetDimension(fem, &N_b));
5655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetNumComponents(fem, &N_comp));
5665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetResidual(prob, field, &f0_func, &f1_func));
5675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetTileSizes(fem, NULL, &N_bl, &N_bc, &N_cb));
56820cf1dd8SToby Isaac   N_bt  = N_b*N_comp;
56920cf1dd8SToby Isaac   N_bst = N_bt*N_q;
57020cf1dd8SToby Isaac   N_t   = N_bst*N_bl;
5715f80ce2aSJacob Faibussowitsch   PetscCheck(N_bc*N_comp == N_t,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of threads %d should be %d * %d", N_t, N_bc, N_comp);
57220cf1dd8SToby Isaac   /* Calculate layout */
57320cf1dd8SToby Isaac   if (Ne % (N_cb*N_bc)) { /* Remainder cells */
5745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEIntegrateResidual_Basic(prob, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
57520cf1dd8SToby Isaac     PetscFunctionReturn(0);
57620cf1dd8SToby Isaac   }
5775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEOpenCLCalculateGrid(fem, Ne, N_cb*N_bc, &x, &y, &z));
57820cf1dd8SToby Isaac   local_work_size[0]  = N_bc*N_comp;
57920cf1dd8SToby Isaac   local_work_size[1]  = 1;
58020cf1dd8SToby Isaac   local_work_size[2]  = 1;
58120cf1dd8SToby Isaac   global_work_size[0] = x * local_work_size[0];
58220cf1dd8SToby Isaac   global_work_size[1] = y * local_work_size[1];
58320cf1dd8SToby Isaac   global_work_size[2] = z * local_work_size[2];
5845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(fem, "GPU layout grid(%d,%d,%d) block(%d,%d,%d) with %d batches\n", x, y, z, local_work_size[0], local_work_size[1], local_work_size[2], N_cb));
5855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(fem, " N_t: %d, N_cb: %d\n", N_t, N_cb));
58620cf1dd8SToby Isaac   /* Generate code */
58720cf1dd8SToby Isaac   if (probAux) {
58820cf1dd8SToby Isaac     PetscSpace P;
58920cf1dd8SToby Isaac     PetscInt   NfAux, order, f;
59020cf1dd8SToby Isaac 
5915f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetNumFields(probAux, &NfAux));
59220cf1dd8SToby Isaac     for (f = 0; f < NfAux; ++f) {
59320cf1dd8SToby Isaac       PetscFE feAux;
59420cf1dd8SToby Isaac 
5955f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSGetDiscretization(probAux, f, (PetscObject *) &feAux));
5965f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFEGetBasisSpace(feAux, &P));
5975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSpaceGetDegree(P, &order, NULL));
5985f80ce2aSJacob Faibussowitsch       PetscCheck(order <= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Can only handle P0 coefficient fields");
59920cf1dd8SToby Isaac     }
60020cf1dd8SToby Isaac   }
6015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEOpenCLGetIntegrationKernel(fem, useAux, &ocl_prog, &ocl_kernel));
60220cf1dd8SToby Isaac   /* Create buffers on the device and send data over */
6035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDataTypeGetSize(ocl->realType, &realSize));
6045f80ce2aSJacob Faibussowitsch   PetscCheck(cgeom->numPoints <= 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support affine geometry for OpenCL integration right now");
60520cf1dd8SToby Isaac   if (sizeof(PetscReal) != realSize) {
60620cf1dd8SToby Isaac     switch (ocl->realType) {
60720cf1dd8SToby Isaac     case PETSC_FLOAT:
60820cf1dd8SToby Isaac     {
60920cf1dd8SToby Isaac       PetscInt c, b, d;
61020cf1dd8SToby Isaac 
6115f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc4(Ne*N_bt,&f_coeff,Ne,&f_coeffAux,Ne*dim*dim,&f_invJ,Ne,&f_detJ));
61220cf1dd8SToby Isaac       for (c = 0; c < Ne; ++c) {
6137be5e748SToby Isaac         f_detJ[c] = (float) cgeom->detJ[c];
61420cf1dd8SToby Isaac         for (d = 0; d < dim*dim; ++d) {
6157be5e748SToby Isaac           f_invJ[c*dim*dim+d] = (float) cgeom->invJ[c * dim * dim + d];
61620cf1dd8SToby Isaac         }
61720cf1dd8SToby Isaac         for (b = 0; b < N_bt; ++b) {
61820cf1dd8SToby Isaac           f_coeff[c*N_bt+b] = (float) coefficients[c*N_bt+b];
61920cf1dd8SToby Isaac         }
62020cf1dd8SToby Isaac       }
62120cf1dd8SToby Isaac       if (coefficientsAux) { /* Assume P0 */
62220cf1dd8SToby Isaac         for (c = 0; c < Ne; ++c) {
62320cf1dd8SToby Isaac           f_coeffAux[c] = (float) coefficientsAux[c];
62420cf1dd8SToby Isaac         }
62520cf1dd8SToby Isaac       }
62620cf1dd8SToby Isaac       oclCoeff      = (void *) f_coeff;
62720cf1dd8SToby Isaac       if (coefficientsAux) {
62820cf1dd8SToby Isaac         oclCoeffAux = (void *) f_coeffAux;
62920cf1dd8SToby Isaac       } else {
63020cf1dd8SToby Isaac         oclCoeffAux = NULL;
63120cf1dd8SToby Isaac       }
63220cf1dd8SToby Isaac       oclInvJ       = (void *) f_invJ;
63320cf1dd8SToby Isaac       oclDetJ       = (void *) f_detJ;
63420cf1dd8SToby Isaac     }
63520cf1dd8SToby Isaac     break;
63620cf1dd8SToby Isaac     case PETSC_DOUBLE:
63720cf1dd8SToby Isaac     {
63820cf1dd8SToby Isaac       PetscInt c, b, d;
63920cf1dd8SToby Isaac 
6405f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc4(Ne*N_bt,&d_coeff,Ne,&d_coeffAux,Ne*dim*dim,&d_invJ,Ne,&d_detJ));
64120cf1dd8SToby Isaac       for (c = 0; c < Ne; ++c) {
6427be5e748SToby Isaac         d_detJ[c] = (double) cgeom->detJ[c];
64320cf1dd8SToby Isaac         for (d = 0; d < dim*dim; ++d) {
6447be5e748SToby Isaac           d_invJ[c*dim*dim+d] = (double) cgeom->invJ[c * dim * dim + d];
64520cf1dd8SToby Isaac         }
64620cf1dd8SToby Isaac         for (b = 0; b < N_bt; ++b) {
64720cf1dd8SToby Isaac           d_coeff[c*N_bt+b] = (double) coefficients[c*N_bt+b];
64820cf1dd8SToby Isaac         }
64920cf1dd8SToby Isaac       }
65020cf1dd8SToby Isaac       if (coefficientsAux) { /* Assume P0 */
65120cf1dd8SToby Isaac         for (c = 0; c < Ne; ++c) {
65220cf1dd8SToby Isaac           d_coeffAux[c] = (double) coefficientsAux[c];
65320cf1dd8SToby Isaac         }
65420cf1dd8SToby Isaac       }
65520cf1dd8SToby Isaac       oclCoeff      = (void *) d_coeff;
65620cf1dd8SToby Isaac       if (coefficientsAux) {
65720cf1dd8SToby Isaac         oclCoeffAux = (void *) d_coeffAux;
65820cf1dd8SToby Isaac       } else {
65920cf1dd8SToby Isaac         oclCoeffAux = NULL;
66020cf1dd8SToby Isaac       }
66120cf1dd8SToby Isaac       oclInvJ       = (void *) d_invJ;
66220cf1dd8SToby Isaac       oclDetJ       = (void *) d_detJ;
66320cf1dd8SToby Isaac     }
66420cf1dd8SToby Isaac     break;
66520cf1dd8SToby Isaac     default:
66698921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported PETSc type %d", ocl->realType);
66720cf1dd8SToby Isaac     }
66820cf1dd8SToby Isaac   } else {
66920cf1dd8SToby Isaac     PetscInt c, d;
67020cf1dd8SToby Isaac 
6715f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc2(Ne*dim*dim,&r_invJ,Ne,&r_detJ));
67220cf1dd8SToby Isaac     for (c = 0; c < Ne; ++c) {
6737be5e748SToby Isaac       r_detJ[c] = cgeom->detJ[c];
67420cf1dd8SToby Isaac       for (d = 0; d < dim*dim; ++d) {
6757be5e748SToby Isaac         r_invJ[c*dim*dim+d] = cgeom->invJ[c * dim * dim + d];
67620cf1dd8SToby Isaac       }
67720cf1dd8SToby Isaac     }
67820cf1dd8SToby Isaac     oclCoeff    = (void *) coefficients;
67920cf1dd8SToby Isaac     oclCoeffAux = (void *) coefficientsAux;
68020cf1dd8SToby Isaac     oclInvJ     = (void *) r_invJ;
68120cf1dd8SToby Isaac     oclDetJ     = (void *) r_detJ;
68220cf1dd8SToby Isaac   }
6832f613bf5SBarry Smith   o_coefficients         = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne*N_bt    * realSize, oclCoeff,    &ierr);
68420cf1dd8SToby Isaac   if (coefficientsAux) {
6852f613bf5SBarry Smith     o_coefficientsAux    = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne         * realSize, oclCoeffAux, &ierr);
68620cf1dd8SToby Isaac   } else {
6872f613bf5SBarry Smith     o_coefficientsAux    = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY,                        Ne         * realSize, oclCoeffAux, &ierr);
68820cf1dd8SToby Isaac   }
6892f613bf5SBarry Smith   o_jacobianInverses     = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne*dim*dim * realSize, oclInvJ,     &ierr);
6902f613bf5SBarry Smith   o_jacobianDeterminants = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne         * realSize, oclDetJ,     &ierr);
6912f613bf5SBarry Smith   o_elemVec              = clCreateBuffer(ocl->ctx_id, CL_MEM_WRITE_ONLY,                       Ne*N_bt    * realSize, NULL,        &ierr);
69220cf1dd8SToby Isaac   /* Kernel launch */
6935f80ce2aSJacob Faibussowitsch   CHKERRQ(clSetKernelArg(ocl_kernel, 0, sizeof(cl_int), (void*) &N_cb));
6945f80ce2aSJacob Faibussowitsch   CHKERRQ(clSetKernelArg(ocl_kernel, 1, sizeof(cl_mem), (void*) &o_coefficients));
6955f80ce2aSJacob Faibussowitsch   CHKERRQ(clSetKernelArg(ocl_kernel, 2, sizeof(cl_mem), (void*) &o_coefficientsAux));
6965f80ce2aSJacob Faibussowitsch   CHKERRQ(clSetKernelArg(ocl_kernel, 3, sizeof(cl_mem), (void*) &o_jacobianInverses));
6975f80ce2aSJacob Faibussowitsch   CHKERRQ(clSetKernelArg(ocl_kernel, 4, sizeof(cl_mem), (void*) &o_jacobianDeterminants));
6985f80ce2aSJacob Faibussowitsch   CHKERRQ(clSetKernelArg(ocl_kernel, 5, sizeof(cl_mem), (void*) &o_elemVec));
6995f80ce2aSJacob Faibussowitsch   CHKERRQ(clEnqueueNDRangeKernel(ocl->queue_id, ocl_kernel, 3, NULL, global_work_size, local_work_size, 0, NULL, &ocl_ev));
70020cf1dd8SToby Isaac   /* Read data back from device */
70120cf1dd8SToby Isaac   if (sizeof(PetscReal) != realSize) {
70220cf1dd8SToby Isaac     switch (ocl->realType) {
70320cf1dd8SToby Isaac     case PETSC_FLOAT:
70420cf1dd8SToby Isaac     {
70520cf1dd8SToby Isaac       float   *elem;
70620cf1dd8SToby Isaac       PetscInt c, b;
70720cf1dd8SToby Isaac 
7085f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree4(f_coeff,f_coeffAux,f_invJ,f_detJ));
7095f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(Ne*N_bt, &elem));
7105f80ce2aSJacob Faibussowitsch       CHKERRQ(clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne*N_bt * realSize, elem, 0, NULL, NULL));
71120cf1dd8SToby Isaac       for (c = 0; c < Ne; ++c) {
71220cf1dd8SToby Isaac         for (b = 0; b < N_bt; ++b) {
71320cf1dd8SToby Isaac           elemVec[c*N_bt+b] = (PetscScalar) elem[c*N_bt+b];
71420cf1dd8SToby Isaac         }
71520cf1dd8SToby Isaac       }
7165f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(elem));
71720cf1dd8SToby Isaac     }
71820cf1dd8SToby Isaac     break;
71920cf1dd8SToby Isaac     case PETSC_DOUBLE:
72020cf1dd8SToby Isaac     {
72120cf1dd8SToby Isaac       double  *elem;
72220cf1dd8SToby Isaac       PetscInt c, b;
72320cf1dd8SToby Isaac 
7245f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree4(d_coeff,d_coeffAux,d_invJ,d_detJ));
7255f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(Ne*N_bt, &elem));
7265f80ce2aSJacob Faibussowitsch       CHKERRQ(clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne*N_bt * realSize, elem, 0, NULL, NULL));
72720cf1dd8SToby Isaac       for (c = 0; c < Ne; ++c) {
72820cf1dd8SToby Isaac         for (b = 0; b < N_bt; ++b) {
72920cf1dd8SToby Isaac           elemVec[c*N_bt+b] = (PetscScalar) elem[c*N_bt+b];
73020cf1dd8SToby Isaac         }
73120cf1dd8SToby Isaac       }
7325f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(elem));
73320cf1dd8SToby Isaac     }
73420cf1dd8SToby Isaac     break;
73520cf1dd8SToby Isaac     default:
73698921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported PETSc type %d", ocl->realType);
73720cf1dd8SToby Isaac     }
73820cf1dd8SToby Isaac   } else {
7395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(r_invJ,r_detJ));
7405f80ce2aSJacob Faibussowitsch     CHKERRQ(clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne*N_bt * realSize, elemVec, 0, NULL, NULL));
74120cf1dd8SToby Isaac   }
74220cf1dd8SToby Isaac   /* Log performance */
7435f80ce2aSJacob Faibussowitsch   CHKERRQ(clGetEventProfilingInfo(ocl_ev, CL_PROFILING_COMMAND_START, sizeof(cl_ulong), &ns_start, NULL));
7445f80ce2aSJacob Faibussowitsch   CHKERRQ(clGetEventProfilingInfo(ocl_ev, CL_PROFILING_COMMAND_END,   sizeof(cl_ulong), &ns_end,   NULL));
74520cf1dd8SToby Isaac   f0Flops = 0;
74620cf1dd8SToby Isaac   switch (ocl->op) {
74720cf1dd8SToby Isaac   case LAPLACIAN:
74820cf1dd8SToby Isaac     f1Flops = useAux ? dim : 0;break;
74920cf1dd8SToby Isaac   case ELASTICITY:
75020cf1dd8SToby Isaac     f1Flops = 2*dim*dim;break;
75120cf1dd8SToby Isaac   }
75220cf1dd8SToby Isaac   numFlops = Ne*(
75320cf1dd8SToby Isaac     N_q*(
75420cf1dd8SToby Isaac       N_b*N_comp*((useField ? 2 : 0) + (useFieldDer ? 2*dim*(dim + 1) : 0))
75520cf1dd8SToby Isaac       /*+
75620cf1dd8SToby Isaac        N_ba*N_compa*((useFieldAux ? 2 : 0) + (useFieldDerAux ? 2*dim*(dim + 1) : 0))*/
75720cf1dd8SToby Isaac       +
75820cf1dd8SToby Isaac       N_comp*((useF0 ? f0Flops + 2 : 0) + (useF1 ? f1Flops + 2*dim : 0)))
75920cf1dd8SToby Isaac     +
76020cf1dd8SToby Isaac     N_b*((useF0 ? 2 : 0) + (useF1 ? 2*dim*(dim + 1) : 0)));
7615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEOpenCLLogResidual(fem, (ns_end - ns_start)*1.0e-9, numFlops));
76220cf1dd8SToby Isaac   /* Cleanup */
7635f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseMemObject(o_coefficients));
7645f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseMemObject(o_coefficientsAux));
7655f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseMemObject(o_jacobianInverses));
7665f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseMemObject(o_jacobianDeterminants));
7675f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseMemObject(o_elemVec));
7685f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseKernel(ocl_kernel));
7695f80ce2aSJacob Faibussowitsch   CHKERRQ(clReleaseProgram(ocl_prog));
77020cf1dd8SToby Isaac   PetscFunctionReturn(0);
77120cf1dd8SToby Isaac }
77220cf1dd8SToby Isaac 
773526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE);
774526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE, PetscInt, const PetscReal [], PetscInt, PetscTabulation);
7752e47ffbbSToby Isaac 
7762b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_OpenCL(PetscFE fem)
77720cf1dd8SToby Isaac {
77820cf1dd8SToby Isaac   PetscFunctionBegin;
77920cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
78020cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
78120cf1dd8SToby Isaac   fem->ops->view                    = NULL;
78220cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_OpenCL;
78320cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
784ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
78520cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_OpenCL;
78620cf1dd8SToby Isaac   fem->ops->integratebdresidual     = NULL/* PetscFEIntegrateBdResidual_OpenCL */;
78720cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_OpenCL */;
78820cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
78920cf1dd8SToby Isaac   PetscFunctionReturn(0);
79020cf1dd8SToby Isaac }
79120cf1dd8SToby Isaac 
79220cf1dd8SToby Isaac /*MC
79320cf1dd8SToby Isaac   PETSCFEOPENCL = "opencl" - A PetscFE object that integrates using a vectorized OpenCL implementation
79420cf1dd8SToby Isaac 
79520cf1dd8SToby Isaac   Level: intermediate
79620cf1dd8SToby Isaac 
79720cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
79820cf1dd8SToby Isaac M*/
79920cf1dd8SToby Isaac 
80020cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_OpenCL(PetscFE fem)
80120cf1dd8SToby Isaac {
80220cf1dd8SToby Isaac   PetscFE_OpenCL *ocl;
80320cf1dd8SToby Isaac   cl_uint         num_platforms;
80420cf1dd8SToby Isaac   cl_platform_id  platform_ids[42];
80520cf1dd8SToby Isaac   cl_uint         num_devices;
80620cf1dd8SToby Isaac   cl_device_id    device_ids[42];
8072da392ccSBarry Smith   cl_int          err;
80820cf1dd8SToby Isaac 
80920cf1dd8SToby Isaac   PetscFunctionBegin;
81020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
8115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(fem,&ocl));
81220cf1dd8SToby Isaac   fem->data = ocl;
81320cf1dd8SToby Isaac 
81420cf1dd8SToby Isaac   /* Init Platform */
8155f80ce2aSJacob Faibussowitsch   CHKERRQ(clGetPlatformIDs(42, platform_ids, &num_platforms));
816*28b400f6SJacob Faibussowitsch   PetscCheck(num_platforms,PetscObjectComm((PetscObject) fem), PETSC_ERR_SUP, "No OpenCL platform found.");
81720cf1dd8SToby Isaac   ocl->pf_id = platform_ids[0];
81820cf1dd8SToby Isaac   /* Init Device */
8195f80ce2aSJacob Faibussowitsch   CHKERRQ(clGetDeviceIDs(ocl->pf_id, CL_DEVICE_TYPE_ALL, 42, device_ids, &num_devices));
820*28b400f6SJacob Faibussowitsch   PetscCheck(num_devices,PetscObjectComm((PetscObject) fem), PETSC_ERR_SUP, "No OpenCL device found.");
82120cf1dd8SToby Isaac   ocl->dev_id = device_ids[0];
82220cf1dd8SToby Isaac   /* Create context with one command queue */
8232da392ccSBarry Smith   ocl->ctx_id   = clCreateContext(0, 1, &(ocl->dev_id), NULL, NULL, &err);CHKERRQ(err);
8242da392ccSBarry Smith   ocl->queue_id = clCreateCommandQueue(ocl->ctx_id, ocl->dev_id, CL_QUEUE_PROFILING_ENABLE, &err);CHKERRQ(err);
82520cf1dd8SToby Isaac   /* Types */
82620cf1dd8SToby Isaac   ocl->realType = PETSC_FLOAT;
82720cf1dd8SToby Isaac   /* Register events */
8285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("OpenCL FEResidual", PETSCFE_CLASSID, &ocl->residualEvent));
82920cf1dd8SToby Isaac   /* Equation handling */
83020cf1dd8SToby Isaac   ocl->op = LAPLACIAN;
83120cf1dd8SToby Isaac 
8325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEInitialize_OpenCL(fem));
83320cf1dd8SToby Isaac   PetscFunctionReturn(0);
83420cf1dd8SToby Isaac }
83520cf1dd8SToby Isaac 
8362b99622eSMatthew G. Knepley /*@
8372b99622eSMatthew G. Knepley   PetscFEOpenCLSetRealType - Set the scalar type for running on the accelerator
8382b99622eSMatthew G. Knepley 
8392b99622eSMatthew G. Knepley   Input Parameters:
8402b99622eSMatthew G. Knepley + fem      - The PetscFE
8412b99622eSMatthew G. Knepley - realType - The scalar type
8422b99622eSMatthew G. Knepley 
8432b99622eSMatthew G. Knepley   Level: developer
8442b99622eSMatthew G. Knepley 
8452b99622eSMatthew G. Knepley .seealso: PetscFEOpenCLGetRealType()
8462b99622eSMatthew G. Knepley @*/
84720cf1dd8SToby Isaac PetscErrorCode PetscFEOpenCLSetRealType(PetscFE fem, PetscDataType realType)
84820cf1dd8SToby Isaac {
84920cf1dd8SToby Isaac   PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
85020cf1dd8SToby Isaac 
85120cf1dd8SToby Isaac   PetscFunctionBegin;
85220cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
85320cf1dd8SToby Isaac   ocl->realType = realType;
85420cf1dd8SToby Isaac   PetscFunctionReturn(0);
85520cf1dd8SToby Isaac }
85620cf1dd8SToby Isaac 
8572b99622eSMatthew G. Knepley /*@
8582b99622eSMatthew G. Knepley   PetscFEOpenCLGetRealType - Get the scalar type for running on the accelerator
8592b99622eSMatthew G. Knepley 
8602b99622eSMatthew G. Knepley   Input Parameter:
8612b99622eSMatthew G. Knepley . fem      - The PetscFE
8622b99622eSMatthew G. Knepley 
8632b99622eSMatthew G. Knepley   Output Parameter:
8642b99622eSMatthew G. Knepley . realType - The scalar type
8652b99622eSMatthew G. Knepley 
8662b99622eSMatthew G. Knepley   Level: developer
8672b99622eSMatthew G. Knepley 
8682b99622eSMatthew G. Knepley .seealso: PetscFEOpenCLSetRealType()
8692b99622eSMatthew G. Knepley @*/
87020cf1dd8SToby Isaac PetscErrorCode PetscFEOpenCLGetRealType(PetscFE fem, PetscDataType *realType)
87120cf1dd8SToby Isaac {
87220cf1dd8SToby Isaac   PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
87320cf1dd8SToby Isaac 
87420cf1dd8SToby Isaac   PetscFunctionBegin;
87520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
87620cf1dd8SToby Isaac   PetscValidPointer(realType, 2);
87720cf1dd8SToby Isaac   *realType = ocl->realType;
87820cf1dd8SToby Isaac   PetscFunctionReturn(0);
87920cf1dd8SToby Isaac }
88020cf1dd8SToby Isaac 
88120cf1dd8SToby Isaac #endif /* PETSC_HAVE_OPENCL */
882