120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac 37be5e748SToby Isaac #if defined(PETSC_HAVE_OPENCL) 420cf1dd8SToby Isaac 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEDestroy_OpenCL(PetscFE fem) 6d71ae5a4SJacob Faibussowitsch { 720cf1dd8SToby Isaac PetscFE_OpenCL *ocl = (PetscFE_OpenCL *)fem->data; 820cf1dd8SToby Isaac 920cf1dd8SToby Isaac PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(clReleaseCommandQueue(ocl->queue_id)); 1120cf1dd8SToby Isaac ocl->queue_id = 0; 129566063dSJacob Faibussowitsch PetscCall(clReleaseContext(ocl->ctx_id)); 1320cf1dd8SToby Isaac ocl->ctx_id = 0; 149566063dSJacob Faibussowitsch PetscCall(PetscFree(ocl)); 1520cf1dd8SToby Isaac PetscFunctionReturn(0); 1620cf1dd8SToby Isaac } 1720cf1dd8SToby Isaac 189371c9d4SSatish Balay #define PetscCallSTR(err) \ 199371c9d4SSatish Balay do { \ 209371c9d4SSatish Balay PetscCall(err); \ 219371c9d4SSatish Balay string_tail += count; \ 229371c9d4SSatish Balay PetscCheck(string_tail != end_of_buffer, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Buffer overflow"); \ 239371c9d4SSatish Balay } while (0) 249371c9d4SSatish Balay enum { 259371c9d4SSatish Balay LAPLACIAN = 0, 269371c9d4SSatish Balay ELASTICITY = 1 279371c9d4SSatish Balay }; 2820cf1dd8SToby Isaac 2920cf1dd8SToby Isaac /* NOTE: This is now broken for vector problems. Must redo loops to respect vector basis elements */ 3020cf1dd8SToby Isaac /* dim Number of spatial dimensions: 2 */ 3120cf1dd8SToby Isaac /* N_b Number of basis functions: generated */ 3220cf1dd8SToby Isaac /* N_{bt} Number of total basis functions: N_b * N_{comp} */ 3320cf1dd8SToby Isaac /* N_q Number of quadrature points: generated */ 3420cf1dd8SToby Isaac /* N_{bs} Number of block cells LCM(N_b, N_q) */ 3520cf1dd8SToby Isaac /* N_{bst} Number of block cell components LCM(N_{bt}, N_q) */ 3620cf1dd8SToby Isaac /* N_{bl} Number of concurrent blocks generated */ 3720cf1dd8SToby Isaac /* N_t Number of threads: N_{bl} * N_{bs} */ 3820cf1dd8SToby Isaac /* N_{cbc} Number of concurrent basis cells: N_{bl} * N_q */ 3920cf1dd8SToby Isaac /* N_{cqc} Number of concurrent quadrature cells: N_{bl} * N_b */ 4020cf1dd8SToby Isaac /* N_{sbc} Number of serial basis cells: N_{bs} / N_q */ 4120cf1dd8SToby Isaac /* N_{sqc} Number of serial quadrature cells: N_{bs} / N_b */ 4220cf1dd8SToby Isaac /* N_{cb} Number of serial cell batches: input */ 4320cf1dd8SToby Isaac /* N_c Number of total cells: N_{cb}*N_{t}/N_{comp} */ 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEOpenCLGenerateIntegrationCode(PetscFE fem, char **string_buffer, PetscInt buffer_length, PetscBool useAux, PetscInt N_bl) 45d71ae5a4SJacob Faibussowitsch { 4620cf1dd8SToby Isaac PetscFE_OpenCL *ocl = (PetscFE_OpenCL *)fem->data; 4720cf1dd8SToby Isaac PetscQuadrature q; 4820cf1dd8SToby Isaac char *string_tail = *string_buffer; 4920cf1dd8SToby Isaac char *end_of_buffer = *string_buffer + buffer_length; 5020cf1dd8SToby Isaac char float_str[] = "float", double_str[] = "double"; 5120cf1dd8SToby Isaac char *numeric_str = &(float_str[0]); 5220cf1dd8SToby Isaac PetscInt op = ocl->op; 5320cf1dd8SToby Isaac PetscBool useField = PETSC_FALSE; 5420cf1dd8SToby Isaac PetscBool useFieldDer = PETSC_TRUE; 5520cf1dd8SToby Isaac PetscBool useFieldAux = useAux; 5620cf1dd8SToby Isaac PetscBool useFieldDerAux = PETSC_FALSE; 5720cf1dd8SToby Isaac PetscBool useF0 = PETSC_TRUE; 5820cf1dd8SToby Isaac PetscBool useF1 = PETSC_TRUE; 5920cf1dd8SToby Isaac const PetscReal *points, *weights; 60ef0bb6c7SMatthew G. Knepley PetscTabulation T; 6120cf1dd8SToby Isaac PetscInt dim, qNc, N_b, N_c, N_q, N_t, p, d, b, c; 6220cf1dd8SToby Isaac size_t count; 6320cf1dd8SToby Isaac 6420cf1dd8SToby Isaac PetscFunctionBegin; 659566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fem, &dim)); 669566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fem, &N_b)); 679566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fem, &N_c)); 689566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fem, &q)); 699566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, &qNc, &N_q, &points, &weights)); 705f80ce2aSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 7120cf1dd8SToby Isaac N_t = N_b * N_c * N_q * N_bl; 7220cf1dd8SToby Isaac /* Enable device extension for double precision */ 7320cf1dd8SToby Isaac if (ocl->realType == PETSC_DOUBLE) { 749566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 7520cf1dd8SToby Isaac "#if defined(cl_khr_fp64)\n" 7620cf1dd8SToby Isaac "# pragma OPENCL EXTENSION cl_khr_fp64: enable\n" 7720cf1dd8SToby Isaac "#elif defined(cl_amd_fp64)\n" 7820cf1dd8SToby Isaac "# pragma OPENCL EXTENSION cl_amd_fp64: enable\n" 7920cf1dd8SToby Isaac "#endif\n", 805f80ce2aSJacob Faibussowitsch &count)); 8120cf1dd8SToby Isaac numeric_str = &(double_str[0]); 8220cf1dd8SToby Isaac } 8320cf1dd8SToby Isaac /* Kernel API */ 849566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 8520cf1dd8SToby Isaac "\n" 8620cf1dd8SToby Isaac "__kernel void integrateElementQuadrature(int N_cb, __global %s *coefficients, __global %s *coefficientsAux, __global %s *jacobianInverses, __global %s *jacobianDeterminants, __global %s *elemVec)\n" 8720cf1dd8SToby Isaac "{\n", 885f80ce2aSJacob Faibussowitsch &count, numeric_str, numeric_str, numeric_str, numeric_str, numeric_str)); 8920cf1dd8SToby Isaac /* Quadrature */ 909566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 9120cf1dd8SToby Isaac " /* Quadrature points\n" 9220cf1dd8SToby Isaac " - (x1,y1,x2,y2,...) */\n" 9320cf1dd8SToby Isaac " const %s points[%d] = {\n", 945f80ce2aSJacob Faibussowitsch &count, numeric_str, N_q * dim)); 9520cf1dd8SToby Isaac for (p = 0; p < N_q; ++p) { 9648a46eb9SPierre Jolivet for (d = 0; d < dim; ++d) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, points[p * dim + d])); 9720cf1dd8SToby Isaac } 989566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count)); 999566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 10020cf1dd8SToby Isaac " /* Quadrature weights\n" 10120cf1dd8SToby Isaac " - (v1,v2,...) */\n" 10220cf1dd8SToby Isaac " const %s weights[%d] = {\n", 1035f80ce2aSJacob Faibussowitsch &count, numeric_str, N_q)); 10448a46eb9SPierre Jolivet for (p = 0; p < N_q; ++p) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, weights[p])); 1059566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count)); 10620cf1dd8SToby Isaac /* Basis Functions */ 1079566063dSJacob Faibussowitsch PetscCall(PetscFEGetCellTabulation(fem, 1, &T)); 1089566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 10920cf1dd8SToby Isaac " /* Nodal basis function evaluations\n" 11020cf1dd8SToby Isaac " - basis component is fastest varying, the basis function, then point */\n" 11120cf1dd8SToby Isaac " const %s Basis[%d] = {\n", 1125f80ce2aSJacob Faibussowitsch &count, numeric_str, N_q * N_b * N_c)); 11320cf1dd8SToby Isaac for (p = 0; p < N_q; ++p) { 11420cf1dd8SToby Isaac for (b = 0; b < N_b; ++b) { 11548a46eb9SPierre Jolivet for (c = 0; c < N_c; ++c) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, T->T[0][(p * N_b + b) * N_c + c])); 11620cf1dd8SToby Isaac } 11720cf1dd8SToby Isaac } 1189566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count)); 1199566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 12020cf1dd8SToby Isaac "\n" 12120cf1dd8SToby Isaac " /* Nodal basis function derivative evaluations,\n" 12220cf1dd8SToby Isaac " - derivative direction is fastest varying, then basis component, then basis function, then point */\n" 12320cf1dd8SToby Isaac " const %s%d BasisDerivatives[%d] = {\n", 1245f80ce2aSJacob Faibussowitsch &count, numeric_str, dim, N_q * N_b * N_c)); 12520cf1dd8SToby Isaac for (p = 0; p < N_q; ++p) { 12620cf1dd8SToby Isaac for (b = 0; b < N_b; ++b) { 12720cf1dd8SToby Isaac for (c = 0; c < N_c; ++c) { 1289566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "(%s%d)(", &count, numeric_str, dim)); 12920cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 13020cf1dd8SToby Isaac if (d > 0) { 1319566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, ", %g", &count, T->T[1][((p * N_b + b) * dim + d) * N_c + c])); 13220cf1dd8SToby Isaac } else { 1339566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g", &count, T->T[1][((p * N_b + b) * dim + d) * N_c + c])); 13420cf1dd8SToby Isaac } 13520cf1dd8SToby Isaac } 1369566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "),\n", &count)); 13720cf1dd8SToby Isaac } 13820cf1dd8SToby Isaac } 13920cf1dd8SToby Isaac } 1409566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count)); 14120cf1dd8SToby Isaac /* Sizes */ 1429566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 14320cf1dd8SToby Isaac " const int dim = %d; // The spatial dimension\n" 14420cf1dd8SToby Isaac " const int N_bl = %d; // The number of concurrent blocks\n" 14520cf1dd8SToby Isaac " const int N_b = %d; // The number of basis functions\n" 14620cf1dd8SToby Isaac " const int N_comp = %d; // The number of basis function components\n" 14720cf1dd8SToby Isaac " const int N_bt = N_b*N_comp; // The total number of scalar basis functions\n" 14820cf1dd8SToby Isaac " const int N_q = %d; // The number of quadrature points\n" 14920cf1dd8SToby 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" 15020cf1dd8SToby Isaac " const int N_t = N_bst*N_bl; // The number of threads, N_bst * N_bl\n" 15120cf1dd8SToby Isaac " const int N_bc = N_t/N_comp; // The number of cells per batch (N_b*N_q*N_bl)\n" 15220cf1dd8SToby Isaac " const int N_sbc = N_bst / (N_q * N_comp);\n" 15320cf1dd8SToby Isaac " const int N_sqc = N_bst / N_bt;\n" 15420cf1dd8SToby Isaac " /*const int N_c = N_cb * N_bc;*/\n" 15520cf1dd8SToby Isaac "\n" 15620cf1dd8SToby Isaac " /* Calculated indices */\n" 15720cf1dd8SToby Isaac " /*const int tidx = get_local_id(0) + get_local_size(0)*get_local_id(1);*/\n" 15820cf1dd8SToby Isaac " const int tidx = get_local_id(0);\n" 15920cf1dd8SToby Isaac " const int blidx = tidx / N_bst; // Block number for this thread\n" 16020cf1dd8SToby Isaac " const int bidx = tidx %% N_bt; // Basis function mapped to this thread\n" 16120cf1dd8SToby Isaac " const int cidx = tidx %% N_comp; // Basis component mapped to this thread\n" 16220cf1dd8SToby Isaac " const int qidx = tidx %% N_q; // Quadrature point mapped to this thread\n" 16320cf1dd8SToby Isaac " const int blbidx = tidx %% N_q + blidx*N_q; // Cell mapped to this thread in the basis phase\n" 16420cf1dd8SToby Isaac " const int blqidx = tidx %% N_b + blidx*N_b; // Cell mapped to this thread in the quadrature phase\n" 16520cf1dd8SToby Isaac " const int gidx = get_group_id(1)*get_num_groups(0) + get_group_id(0);\n" 16620cf1dd8SToby Isaac " const int Goffset = gidx*N_cb*N_bc;\n", 1675f80ce2aSJacob Faibussowitsch &count, dim, N_bl, N_b, N_c, N_q)); 16820cf1dd8SToby Isaac /* Local memory */ 1699566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 17020cf1dd8SToby Isaac "\n" 17120cf1dd8SToby Isaac " /* Quadrature data */\n" 17220cf1dd8SToby Isaac " %s w; // $w_q$, Quadrature weight at $x_q$\n" 17320cf1dd8SToby Isaac " __local %s phi_i[%d]; //[N_bt*N_q]; // $\\phi_i(x_q)$, Value of the basis function $i$ at $x_q$\n" 17420cf1dd8SToby 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" 17520cf1dd8SToby Isaac " /* Geometric data */\n" 17620cf1dd8SToby Isaac " __local %s detJ[%d]; //[N_t]; // $|J(x_q)|$, Jacobian determinant at $x_q$\n" 17720cf1dd8SToby Isaac " __local %s invJ[%d];//[N_t*dim*dim]; // $J^{-1}(x_q)$, Jacobian inverse at $x_q$\n", 1789371c9d4SSatish Balay &count, numeric_str, numeric_str, N_b * N_c * N_q, numeric_str, dim, N_b * N_c * N_q, numeric_str, N_t, numeric_str, N_t * dim * dim)); 1799566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 18020cf1dd8SToby Isaac " /* FEM data */\n" 18120cf1dd8SToby 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", 1825f80ce2aSJacob Faibussowitsch &count, numeric_str, N_t * N_b * N_c)); 18320cf1dd8SToby Isaac if (useAux) { 1849371c9d4SSatish Balay PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " __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", &count, numeric_str, N_t)); 18520cf1dd8SToby Isaac } 18620cf1dd8SToby Isaac if (useF0) { 1879566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 18820cf1dd8SToby Isaac " /* Intermediate calculations */\n" 18920cf1dd8SToby 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", 1905f80ce2aSJacob Faibussowitsch &count, numeric_str, N_t * N_q)); 19120cf1dd8SToby Isaac } 19248a46eb9SPierre Jolivet if (useF1) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " __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", &count, numeric_str, dim, N_t * N_q)); 19320cf1dd8SToby Isaac /* TODO: If using elasticity, put in mu/lambda coefficients */ 1949566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 19520cf1dd8SToby Isaac " /* Output data */\n" 19620cf1dd8SToby Isaac " %s e_i; // Coefficient $e_i$ of the residual\n\n", 1975f80ce2aSJacob Faibussowitsch &count, numeric_str)); 19820cf1dd8SToby Isaac /* One-time loads */ 1999566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 20020cf1dd8SToby Isaac " /* These should be generated inline */\n" 20120cf1dd8SToby Isaac " /* Load quadrature weights */\n" 20220cf1dd8SToby Isaac " w = weights[qidx];\n" 20320cf1dd8SToby Isaac " /* Load basis tabulation \\phi_i for this cell */\n" 20420cf1dd8SToby Isaac " if (tidx < N_bt*N_q) {\n" 20520cf1dd8SToby Isaac " phi_i[tidx] = Basis[tidx];\n" 20620cf1dd8SToby Isaac " phiDer_i[tidx] = BasisDerivatives[tidx];\n" 20720cf1dd8SToby Isaac " }\n\n", 2085f80ce2aSJacob Faibussowitsch &count)); 20920cf1dd8SToby Isaac /* Batch loads */ 2109566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 21120cf1dd8SToby Isaac " for (int batch = 0; batch < N_cb; ++batch) {\n" 21220cf1dd8SToby Isaac " /* Load geometry */\n" 21320cf1dd8SToby Isaac " detJ[tidx] = jacobianDeterminants[Goffset+batch*N_bc+tidx];\n" 21420cf1dd8SToby Isaac " for (int n = 0; n < dim*dim; ++n) {\n" 21520cf1dd8SToby Isaac " const int offset = n*N_t;\n" 21620cf1dd8SToby Isaac " invJ[offset+tidx] = jacobianInverses[(Goffset+batch*N_bc)*dim*dim+offset+tidx];\n" 21720cf1dd8SToby Isaac " }\n" 21820cf1dd8SToby Isaac " /* Load coefficients u_i for this cell */\n" 21920cf1dd8SToby Isaac " for (int n = 0; n < N_bt; ++n) {\n" 22020cf1dd8SToby Isaac " const int offset = n*N_t;\n" 22120cf1dd8SToby Isaac " u_i[offset+tidx] = coefficients[(Goffset*N_bt)+batch*N_t*N_b+offset+tidx];\n" 22220cf1dd8SToby Isaac " }\n", 2235f80ce2aSJacob Faibussowitsch &count)); 22420cf1dd8SToby Isaac if (useAux) { 2259566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 22620cf1dd8SToby Isaac " /* Load coefficients a_i for this cell */\n" 22720cf1dd8SToby Isaac " /* TODO: This should not be N_t here, it should be N_bc*N_comp_aux */\n" 22820cf1dd8SToby Isaac " a_i[tidx] = coefficientsAux[Goffset+batch*N_t+tidx];\n", 2295f80ce2aSJacob Faibussowitsch &count)); 23020cf1dd8SToby Isaac } 23120cf1dd8SToby Isaac /* Quadrature phase */ 2329566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 23320cf1dd8SToby Isaac " barrier(CLK_LOCAL_MEM_FENCE);\n" 23420cf1dd8SToby Isaac "\n" 23520cf1dd8SToby Isaac " /* Map coefficients to values at quadrature points */\n" 23620cf1dd8SToby Isaac " for (int c = 0; c < N_sqc; ++c) {\n" 23720cf1dd8SToby Isaac " const int cell = c*N_bl*N_b + blqidx;\n" 23820cf1dd8SToby Isaac " const int fidx = (cell*N_q + qidx)*N_comp + cidx;\n", 2395f80ce2aSJacob Faibussowitsch &count)); 24048a46eb9SPierre Jolivet if (useField) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " %s u[%d]; //[N_comp]; // $u(x_q)$, Value of the field at $x_q$\n", &count, numeric_str, N_c)); 24148a46eb9SPierre Jolivet if (useFieldDer) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " %s%d gradU[%d]; //[N_comp]; // $\\nabla u(x_q)$, Value of the field gradient at $x_q$\n", &count, numeric_str, dim, N_c)); 24248a46eb9SPierre Jolivet if (useFieldAux) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " %s a[%d]; //[1]; // $a(x_q)$, Value of the auxiliary fields at $x_q$\n", &count, numeric_str, 1)); 24348a46eb9SPierre Jolivet if (useFieldDerAux) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " %s%d gradA[%d]; //[1]; // $\\nabla a(x_q)$, Value of the auxiliary field gradient at $x_q$\n", &count, numeric_str, dim, 1)); 2449566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 24520cf1dd8SToby Isaac "\n" 24620cf1dd8SToby Isaac " for (int comp = 0; comp < N_comp; ++comp) {\n", 2475f80ce2aSJacob Faibussowitsch &count)); 2489566063dSJacob Faibussowitsch if (useField) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " u[comp] = 0.0;\n", &count)); 24920cf1dd8SToby Isaac if (useFieldDer) { 25020cf1dd8SToby Isaac switch (dim) { 251d71ae5a4SJacob Faibussowitsch case 1: 252d71ae5a4SJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradU[comp].x = 0.0;\n", &count)); 253d71ae5a4SJacob Faibussowitsch break; 254d71ae5a4SJacob Faibussowitsch case 2: 255d71ae5a4SJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradU[comp].x = 0.0; gradU[comp].y = 0.0;\n", &count)); 256d71ae5a4SJacob Faibussowitsch break; 257d71ae5a4SJacob Faibussowitsch case 3: 258d71ae5a4SJacob Faibussowitsch PetscCallSTR(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)); 259d71ae5a4SJacob Faibussowitsch break; 26020cf1dd8SToby Isaac } 26120cf1dd8SToby Isaac } 2629371c9d4SSatish Balay PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " }\n", &count)); 26348a46eb9SPierre Jolivet if (useFieldAux) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " a[0] = 0.0;\n", &count)); 26420cf1dd8SToby Isaac if (useFieldDerAux) { 26520cf1dd8SToby Isaac switch (dim) { 266d71ae5a4SJacob Faibussowitsch case 1: 267d71ae5a4SJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradA[0].x = 0.0;\n", &count)); 268d71ae5a4SJacob Faibussowitsch break; 269d71ae5a4SJacob Faibussowitsch case 2: 270d71ae5a4SJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradA[0].x = 0.0; gradA[0].y = 0.0;\n", &count)); 271d71ae5a4SJacob Faibussowitsch break; 272d71ae5a4SJacob Faibussowitsch case 3: 273d71ae5a4SJacob Faibussowitsch PetscCallSTR(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)); 274d71ae5a4SJacob Faibussowitsch break; 27520cf1dd8SToby Isaac } 27620cf1dd8SToby Isaac } 2779566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 27820cf1dd8SToby Isaac " /* Get field and derivatives at this quadrature point */\n" 27920cf1dd8SToby Isaac " for (int i = 0; i < N_b; ++i) {\n" 28020cf1dd8SToby Isaac " for (int comp = 0; comp < N_comp; ++comp) {\n" 28120cf1dd8SToby Isaac " const int b = i*N_comp+comp;\n" 28220cf1dd8SToby Isaac " const int pidx = qidx*N_bt + b;\n" 28320cf1dd8SToby Isaac " const int uidx = cell*N_bt + b;\n" 28420cf1dd8SToby Isaac " %s%d realSpaceDer;\n\n", 2855f80ce2aSJacob Faibussowitsch &count, numeric_str, dim)); 2869566063dSJacob Faibussowitsch if (useField) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " u[comp] += u_i[uidx]*phi_i[pidx];\n", &count)); 28720cf1dd8SToby Isaac if (useFieldDer) { 28820cf1dd8SToby Isaac switch (dim) { 28920cf1dd8SToby Isaac case 2: 2909566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 29120cf1dd8SToby 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" 29220cf1dd8SToby Isaac " gradU[comp].x += u_i[uidx]*realSpaceDer.x;\n" 29320cf1dd8SToby 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" 29420cf1dd8SToby Isaac " gradU[comp].y += u_i[uidx]*realSpaceDer.y;\n", 2959371c9d4SSatish Balay &count)); 2969371c9d4SSatish Balay break; 29720cf1dd8SToby Isaac case 3: 2989566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 29920cf1dd8SToby 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" 30020cf1dd8SToby Isaac " gradU[comp].x += u_i[uidx]*realSpaceDer.x;\n" 30120cf1dd8SToby 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" 30220cf1dd8SToby Isaac " gradU[comp].y += u_i[uidx]*realSpaceDer.y;\n" 30320cf1dd8SToby 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" 30420cf1dd8SToby Isaac " gradU[comp].z += u_i[uidx]*realSpaceDer.z;\n", 3059371c9d4SSatish Balay &count)); 3069371c9d4SSatish Balay break; 30720cf1dd8SToby Isaac } 30820cf1dd8SToby Isaac } 3099566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 31020cf1dd8SToby Isaac " }\n" 31120cf1dd8SToby Isaac " }\n", 3125f80ce2aSJacob Faibussowitsch &count)); 31348a46eb9SPierre Jolivet if (useFieldAux) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " a[0] += a_i[cell];\n", &count)); 31420cf1dd8SToby Isaac /* Calculate residual at quadrature points: Should be generated by an weak form egine */ 3159371c9d4SSatish Balay PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " /* Process values at quadrature points */\n", &count)); 31620cf1dd8SToby Isaac switch (op) { 31720cf1dd8SToby Isaac case LAPLACIAN: 31848a46eb9SPierre Jolivet if (useF0) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_0[fidx] = 4.0;\n", &count)); 31920cf1dd8SToby Isaac if (useF1) { 3209566063dSJacob Faibussowitsch if (useAux) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_1[fidx] = a[0]*gradU[cidx];\n", &count)); 3219566063dSJacob Faibussowitsch else PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_1[fidx] = gradU[cidx];\n", &count)); 32220cf1dd8SToby Isaac } 32320cf1dd8SToby Isaac break; 32420cf1dd8SToby Isaac case ELASTICITY: 3259566063dSJacob Faibussowitsch if (useF0) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_0[fidx] = 4.0;\n", &count)); 32620cf1dd8SToby Isaac if (useF1) { 32720cf1dd8SToby Isaac switch (dim) { 32820cf1dd8SToby Isaac case 2: 3299566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 33020cf1dd8SToby Isaac " switch (cidx) {\n" 33120cf1dd8SToby Isaac " case 0:\n" 33220cf1dd8SToby Isaac " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[0].x + gradU[0].x);\n" 33320cf1dd8SToby Isaac " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[0].y + gradU[1].x);\n" 33420cf1dd8SToby Isaac " break;\n" 33520cf1dd8SToby Isaac " case 1:\n" 33620cf1dd8SToby Isaac " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[1].x + gradU[0].y);\n" 33720cf1dd8SToby Isaac " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[1].y + gradU[1].y);\n" 33820cf1dd8SToby Isaac " }\n", 3399371c9d4SSatish Balay &count)); 3409371c9d4SSatish Balay break; 34120cf1dd8SToby Isaac case 3: 3429566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 34320cf1dd8SToby Isaac " switch (cidx) {\n" 34420cf1dd8SToby Isaac " case 0:\n" 34520cf1dd8SToby Isaac " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].x + gradU[0].x);\n" 34620cf1dd8SToby Isaac " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].y + gradU[1].x);\n" 34720cf1dd8SToby Isaac " f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].z + gradU[2].x);\n" 34820cf1dd8SToby Isaac " break;\n" 34920cf1dd8SToby Isaac " case 1:\n" 35020cf1dd8SToby Isaac " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].x + gradU[0].y);\n" 35120cf1dd8SToby Isaac " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].y + gradU[1].y);\n" 35220cf1dd8SToby Isaac " f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].y + gradU[2].y);\n" 35320cf1dd8SToby Isaac " break;\n" 35420cf1dd8SToby Isaac " case 2:\n" 35520cf1dd8SToby Isaac " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].x + gradU[0].z);\n" 35620cf1dd8SToby Isaac " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].y + gradU[1].z);\n" 35720cf1dd8SToby Isaac " f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].y + gradU[2].z);\n" 35820cf1dd8SToby Isaac " }\n", 3599371c9d4SSatish Balay &count)); 36020cf1dd8SToby Isaac break; 3619371c9d4SSatish Balay } 3629371c9d4SSatish Balay } 3639371c9d4SSatish Balay break; 364d71ae5a4SJacob Faibussowitsch default: 365d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "PDE operator %d is not supported", op); 36620cf1dd8SToby Isaac } 3679566063dSJacob Faibussowitsch if (useF0) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_0[fidx] *= detJ[cell]*w;\n", &count)); 36820cf1dd8SToby Isaac if (useF1) { 36920cf1dd8SToby Isaac switch (dim) { 370d71ae5a4SJacob Faibussowitsch case 1: 371d71ae5a4SJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_1[fidx].x *= detJ[cell]*w;\n", &count)); 372d71ae5a4SJacob Faibussowitsch break; 373d71ae5a4SJacob Faibussowitsch case 2: 374d71ae5a4SJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_1[fidx].x *= detJ[cell]*w; f_1[fidx].y *= detJ[cell]*w;\n", &count)); 375d71ae5a4SJacob Faibussowitsch break; 376d71ae5a4SJacob Faibussowitsch case 3: 377d71ae5a4SJacob Faibussowitsch PetscCallSTR(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)); 378d71ae5a4SJacob Faibussowitsch break; 37920cf1dd8SToby Isaac } 38020cf1dd8SToby Isaac } 38120cf1dd8SToby Isaac /* Thread transpose */ 3829566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 38320cf1dd8SToby Isaac " }\n\n" 38420cf1dd8SToby Isaac " /* ==== TRANSPOSE THREADS ==== */\n" 38520cf1dd8SToby Isaac " barrier(CLK_LOCAL_MEM_FENCE);\n\n", 3865f80ce2aSJacob Faibussowitsch &count)); 38720cf1dd8SToby Isaac /* Basis phase */ 3889566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 38920cf1dd8SToby Isaac " /* Map values at quadrature points to coefficients */\n" 39020cf1dd8SToby Isaac " for (int c = 0; c < N_sbc; ++c) {\n" 39120cf1dd8SToby Isaac " const int cell = c*N_bl*N_q + blbidx; /* Cell number in batch */\n" 39220cf1dd8SToby Isaac "\n" 39320cf1dd8SToby Isaac " e_i = 0.0;\n" 39420cf1dd8SToby Isaac " for (int q = 0; q < N_q; ++q) {\n" 39520cf1dd8SToby Isaac " const int pidx = q*N_bt + bidx;\n" 39620cf1dd8SToby Isaac " const int fidx = (cell*N_q + q)*N_comp + cidx;\n" 39720cf1dd8SToby Isaac " %s%d realSpaceDer;\n\n", 3985f80ce2aSJacob Faibussowitsch &count, numeric_str, dim)); 39920cf1dd8SToby Isaac 4009566063dSJacob Faibussowitsch if (useF0) PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " e_i += phi_i[pidx]*f_0[fidx];\n", &count)); 40120cf1dd8SToby Isaac if (useF1) { 40220cf1dd8SToby Isaac switch (dim) { 40320cf1dd8SToby Isaac case 2: 4049566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 40520cf1dd8SToby 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" 40620cf1dd8SToby Isaac " e_i += realSpaceDer.x*f_1[fidx].x;\n" 40720cf1dd8SToby 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" 40820cf1dd8SToby Isaac " e_i += realSpaceDer.y*f_1[fidx].y;\n", 4099371c9d4SSatish Balay &count)); 4109371c9d4SSatish Balay break; 41120cf1dd8SToby Isaac case 3: 4129566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 41320cf1dd8SToby 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" 41420cf1dd8SToby Isaac " e_i += realSpaceDer.x*f_1[fidx].x;\n" 41520cf1dd8SToby 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" 41620cf1dd8SToby Isaac " e_i += realSpaceDer.y*f_1[fidx].y;\n" 41720cf1dd8SToby 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" 41820cf1dd8SToby Isaac " e_i += realSpaceDer.z*f_1[fidx].z;\n", 4199371c9d4SSatish Balay &count)); 4209371c9d4SSatish Balay break; 42120cf1dd8SToby Isaac } 42220cf1dd8SToby Isaac } 4239566063dSJacob Faibussowitsch PetscCallSTR(PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, 42420cf1dd8SToby Isaac " }\n" 42520cf1dd8SToby Isaac " /* Write element vector for N_{cbc} cells at a time */\n" 42620cf1dd8SToby Isaac " elemVec[(Goffset + batch*N_bc + c*N_bl*N_q)*N_bt + tidx] = e_i;\n" 42720cf1dd8SToby Isaac " }\n" 42820cf1dd8SToby Isaac " /* ==== Could do one write per batch ==== */\n" 42920cf1dd8SToby Isaac " }\n" 43020cf1dd8SToby Isaac " return;\n" 43120cf1dd8SToby Isaac "}\n", 4325f80ce2aSJacob Faibussowitsch &count)); 43320cf1dd8SToby Isaac PetscFunctionReturn(0); 43420cf1dd8SToby Isaac } 43520cf1dd8SToby Isaac 436d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEOpenCLGetIntegrationKernel(PetscFE fem, PetscBool useAux, cl_program *ocl_prog, cl_kernel *ocl_kernel) 437d71ae5a4SJacob Faibussowitsch { 43820cf1dd8SToby Isaac PetscFE_OpenCL *ocl = (PetscFE_OpenCL *)fem->data; 43920cf1dd8SToby Isaac PetscInt dim, N_bl; 44020cf1dd8SToby Isaac PetscBool flg; 44120cf1dd8SToby Isaac char *buffer; 44220cf1dd8SToby Isaac size_t len; 44320cf1dd8SToby Isaac char errMsg[8192]; 4442da392ccSBarry Smith cl_int err; 44520cf1dd8SToby Isaac 44620cf1dd8SToby Isaac PetscFunctionBegin; 4479566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fem, &dim)); 4489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(8192, &buffer)); 4499566063dSJacob Faibussowitsch PetscCall(PetscFEGetTileSizes(fem, NULL, &N_bl, NULL, NULL)); 4509566063dSJacob Faibussowitsch PetscCall(PetscFEOpenCLGenerateIntegrationCode(fem, &buffer, 8192, useAux, N_bl)); 4519566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)fem)->options, ((PetscObject)fem)->prefix, "-petscfe_opencl_kernel_print", &flg)); 4529566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)fem), "OpenCL FE Integration Kernel:\n%s\n", buffer)); 4539566063dSJacob Faibussowitsch PetscCall(PetscStrlen(buffer, &len)); 4549371c9d4SSatish Balay *ocl_prog = clCreateProgramWithSource(ocl->ctx_id, 1, (const char **)&buffer, &len, &err); 4559371c9d4SSatish Balay PetscCall(err); 4562da392ccSBarry Smith err = clBuildProgram(*ocl_prog, 0, NULL, NULL, NULL, NULL); 4572da392ccSBarry Smith if (err != CL_SUCCESS) { 4582f613bf5SBarry Smith err = clGetProgramBuildInfo(*ocl_prog, ocl->dev_id, CL_PROGRAM_BUILD_LOG, 8192 * sizeof(char), &errMsg, NULL); 45998921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Build failed! Log:\n %s", errMsg); 46020cf1dd8SToby Isaac } 4619566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 462d0609cedSBarry Smith *ocl_kernel = clCreateKernel(*ocl_prog, "integrateElementQuadrature", &err); 46320cf1dd8SToby Isaac PetscFunctionReturn(0); 46420cf1dd8SToby Isaac } 46520cf1dd8SToby Isaac 466d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEOpenCLCalculateGrid(PetscFE fem, PetscInt N, PetscInt blockSize, size_t *x, size_t *y, size_t *z) 467d71ae5a4SJacob Faibussowitsch { 46820cf1dd8SToby Isaac const PetscInt Nblocks = N / blockSize; 46920cf1dd8SToby Isaac 47020cf1dd8SToby Isaac PetscFunctionBegin; 4715f80ce2aSJacob Faibussowitsch PetscCheck(!(N % blockSize), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid block size %d for %d elements", blockSize, N); 47220cf1dd8SToby Isaac *z = 1; 473623c9f23SSatish Balay *y = 1; 47420cf1dd8SToby Isaac for (*x = (size_t)(PetscSqrtReal(Nblocks) + 0.5); *x > 0; --*x) { 47520cf1dd8SToby Isaac *y = Nblocks / *x; 476526996e7SStefano Zampini if (*x * *y == (size_t)Nblocks) break; 47720cf1dd8SToby Isaac } 47863a3b9bcSJacob Faibussowitsch PetscCheck(*x * *y == (size_t)Nblocks, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Could not find partition for %" PetscInt_FMT " with block size %" PetscInt_FMT, N, blockSize); 47920cf1dd8SToby Isaac PetscFunctionReturn(0); 48020cf1dd8SToby Isaac } 48120cf1dd8SToby Isaac 482d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEOpenCLLogResidual(PetscFE fem, PetscLogDouble time, PetscLogDouble flops) 483d71ae5a4SJacob Faibussowitsch { 48420cf1dd8SToby Isaac PetscFE_OpenCL *ocl = (PetscFE_OpenCL *)fem->data; 48520cf1dd8SToby Isaac PetscStageLog stageLog; 48620cf1dd8SToby Isaac PetscEventPerfLog eventLog = NULL; 487afbfd3a7SStefano Zampini int stage; 48820cf1dd8SToby Isaac 48920cf1dd8SToby Isaac PetscFunctionBegin; 4909566063dSJacob Faibussowitsch PetscCall(PetscLogGetStageLog(&stageLog)); 4919566063dSJacob Faibussowitsch PetscCall(PetscStageLogGetCurrent(stageLog, &stage)); 4929566063dSJacob Faibussowitsch PetscCall(PetscStageLogGetEventPerfLog(stageLog, stage, &eventLog)); 49320cf1dd8SToby Isaac /* Log performance info */ 49420cf1dd8SToby Isaac eventLog->eventInfo[ocl->residualEvent].count++; 49520cf1dd8SToby Isaac eventLog->eventInfo[ocl->residualEvent].time += time; 49620cf1dd8SToby Isaac eventLog->eventInfo[ocl->residualEvent].flops += flops; 49720cf1dd8SToby Isaac PetscFunctionReturn(0); 49820cf1dd8SToby Isaac } 49920cf1dd8SToby Isaac 500d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEIntegrateResidual_OpenCL(PetscDS prob, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 501d71ae5a4SJacob Faibussowitsch { 50220cf1dd8SToby Isaac /* Nbc = batchSize */ 503849189efSMatthew G. Knepley PetscFE fem; 504849189efSMatthew G. Knepley PetscFE_OpenCL *ocl; 50520cf1dd8SToby Isaac PetscPointFunc f0_func; 50620cf1dd8SToby Isaac PetscPointFunc f1_func; 50720cf1dd8SToby Isaac PetscQuadrature q; 50820cf1dd8SToby Isaac PetscInt dim, qNc; 50920cf1dd8SToby Isaac PetscInt N_b; /* The number of basis functions */ 51020cf1dd8SToby Isaac PetscInt N_comp; /* The number of basis function components */ 51120cf1dd8SToby Isaac PetscInt N_bt; /* The total number of scalar basis functions */ 51220cf1dd8SToby Isaac PetscInt N_q; /* The number of quadrature points */ 51320cf1dd8SToby Isaac PetscInt N_bst; /* The block size, LCM(N_bt, N_q), Notice that a block is not process simultaneously */ 51420cf1dd8SToby Isaac PetscInt N_t; /* The number of threads, N_bst * N_bl */ 51520cf1dd8SToby Isaac PetscInt N_bl; /* The number of blocks */ 51620cf1dd8SToby Isaac PetscInt N_bc; /* The batch size, N_bl*N_q*N_b */ 51720cf1dd8SToby Isaac PetscInt N_cb; /* The number of batches */ 5186528b96dSMatthew G. Knepley const PetscInt field = key.field; 51920cf1dd8SToby Isaac PetscInt numFlops, f0Flops = 0, f1Flops = 0; 52020cf1dd8SToby Isaac PetscBool useAux = probAux ? PETSC_TRUE : PETSC_FALSE; 52120cf1dd8SToby Isaac PetscBool useField = PETSC_FALSE; 52220cf1dd8SToby Isaac PetscBool useFieldDer = PETSC_TRUE; 52320cf1dd8SToby Isaac PetscBool useF0 = PETSC_TRUE; 52420cf1dd8SToby Isaac PetscBool useF1 = PETSC_TRUE; 52520cf1dd8SToby Isaac /* OpenCL variables */ 52620cf1dd8SToby Isaac cl_program ocl_prog; 52720cf1dd8SToby Isaac cl_kernel ocl_kernel; 52820cf1dd8SToby Isaac cl_event ocl_ev; /* The event for tracking kernel execution */ 52920cf1dd8SToby Isaac cl_ulong ns_start; /* Nanoseconds counter on GPU at kernel start */ 53020cf1dd8SToby Isaac cl_ulong ns_end; /* Nanoseconds counter on GPU at kernel stop */ 53120cf1dd8SToby Isaac cl_mem o_jacobianInverses, o_jacobianDeterminants; 53220cf1dd8SToby Isaac cl_mem o_coefficients, o_coefficientsAux, o_elemVec; 53320cf1dd8SToby Isaac float *f_coeff = NULL, *f_coeffAux = NULL, *f_invJ = NULL, *f_detJ = NULL; 53420cf1dd8SToby Isaac double *d_coeff = NULL, *d_coeffAux = NULL, *d_invJ = NULL, *d_detJ = NULL; 53520cf1dd8SToby Isaac PetscReal *r_invJ = NULL, *r_detJ = NULL; 53620cf1dd8SToby Isaac void *oclCoeff, *oclCoeffAux, *oclInvJ, *oclDetJ; 53720cf1dd8SToby Isaac size_t local_work_size[3], global_work_size[3]; 53820cf1dd8SToby Isaac size_t realSize, x, y, z; 53920cf1dd8SToby Isaac const PetscReal *points, *weights; 540d0609cedSBarry Smith int err; 54120cf1dd8SToby Isaac 54220cf1dd8SToby Isaac PetscFunctionBegin; 5439566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fem)); 544849189efSMatthew G. Knepley ocl = (PetscFE_OpenCL *)fem->data; 5459371c9d4SSatish Balay if (!Ne) { 5469371c9d4SSatish Balay PetscCall(PetscFEOpenCLLogResidual(fem, 0.0, 0.0)); 5479371c9d4SSatish Balay PetscFunctionReturn(0); 5489371c9d4SSatish Balay } 5499566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fem, &dim)); 5509566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fem, &q)); 5519566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, &qNc, &N_q, &points, &weights)); 55263a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 5539566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fem, &N_b)); 5549566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fem, &N_comp)); 5559566063dSJacob Faibussowitsch PetscCall(PetscDSGetResidual(prob, field, &f0_func, &f1_func)); 5569566063dSJacob Faibussowitsch PetscCall(PetscFEGetTileSizes(fem, NULL, &N_bl, &N_bc, &N_cb)); 55720cf1dd8SToby Isaac N_bt = N_b * N_comp; 55820cf1dd8SToby Isaac N_bst = N_bt * N_q; 55920cf1dd8SToby Isaac N_t = N_bst * N_bl; 5605f80ce2aSJacob 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); 56120cf1dd8SToby Isaac /* Calculate layout */ 56220cf1dd8SToby Isaac if (Ne % (N_cb * N_bc)) { /* Remainder cells */ 5639566063dSJacob Faibussowitsch PetscCall(PetscFEIntegrateResidual_Basic(prob, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 56420cf1dd8SToby Isaac PetscFunctionReturn(0); 56520cf1dd8SToby Isaac } 5669566063dSJacob Faibussowitsch PetscCall(PetscFEOpenCLCalculateGrid(fem, Ne, N_cb * N_bc, &x, &y, &z)); 56720cf1dd8SToby Isaac local_work_size[0] = N_bc * N_comp; 56820cf1dd8SToby Isaac local_work_size[1] = 1; 56920cf1dd8SToby Isaac local_work_size[2] = 1; 57020cf1dd8SToby Isaac global_work_size[0] = x * local_work_size[0]; 57120cf1dd8SToby Isaac global_work_size[1] = y * local_work_size[1]; 57220cf1dd8SToby Isaac global_work_size[2] = z * local_work_size[2]; 57363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(fem, "GPU layout grid(%zu,%zu,%zu) block(%zu,%zu,%zu) with %d batches\n", x, y, z, local_work_size[0], local_work_size[1], local_work_size[2], N_cb)); 5749566063dSJacob Faibussowitsch PetscCall(PetscInfo(fem, " N_t: %d, N_cb: %d\n", N_t, N_cb)); 57520cf1dd8SToby Isaac /* Generate code */ 57620cf1dd8SToby Isaac if (probAux) { 57720cf1dd8SToby Isaac PetscSpace P; 57820cf1dd8SToby Isaac PetscInt NfAux, order, f; 57920cf1dd8SToby Isaac 5809566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(probAux, &NfAux)); 58120cf1dd8SToby Isaac for (f = 0; f < NfAux; ++f) { 58220cf1dd8SToby Isaac PetscFE feAux; 58320cf1dd8SToby Isaac 5849566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(probAux, f, (PetscObject *)&feAux)); 5859566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(feAux, &P)); 5869566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(P, &order, NULL)); 5875f80ce2aSJacob Faibussowitsch PetscCheck(order <= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Can only handle P0 coefficient fields"); 58820cf1dd8SToby Isaac } 58920cf1dd8SToby Isaac } 5909566063dSJacob Faibussowitsch PetscCall(PetscFEOpenCLGetIntegrationKernel(fem, useAux, &ocl_prog, &ocl_kernel)); 59120cf1dd8SToby Isaac /* Create buffers on the device and send data over */ 5929566063dSJacob Faibussowitsch PetscCall(PetscDataTypeGetSize(ocl->realType, &realSize)); 5935f80ce2aSJacob Faibussowitsch PetscCheck(cgeom->numPoints <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support affine geometry for OpenCL integration right now"); 59420cf1dd8SToby Isaac if (sizeof(PetscReal) != realSize) { 59520cf1dd8SToby Isaac switch (ocl->realType) { 5969371c9d4SSatish Balay case PETSC_FLOAT: { 59720cf1dd8SToby Isaac PetscInt c, b, d; 59820cf1dd8SToby Isaac 5999566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(Ne * N_bt, &f_coeff, Ne, &f_coeffAux, Ne * dim * dim, &f_invJ, Ne, &f_detJ)); 60020cf1dd8SToby Isaac for (c = 0; c < Ne; ++c) { 6017be5e748SToby Isaac f_detJ[c] = (float)cgeom->detJ[c]; 602ad540459SPierre Jolivet for (d = 0; d < dim * dim; ++d) f_invJ[c * dim * dim + d] = (float)cgeom->invJ[c * dim * dim + d]; 603ad540459SPierre Jolivet for (b = 0; b < N_bt; ++b) f_coeff[c * N_bt + b] = (float)coefficients[c * N_bt + b]; 60420cf1dd8SToby Isaac } 60520cf1dd8SToby Isaac if (coefficientsAux) { /* Assume P0 */ 606ad540459SPierre Jolivet for (c = 0; c < Ne; ++c) f_coeffAux[c] = (float)coefficientsAux[c]; 60720cf1dd8SToby Isaac } 60820cf1dd8SToby Isaac oclCoeff = (void *)f_coeff; 60920cf1dd8SToby Isaac if (coefficientsAux) { 61020cf1dd8SToby Isaac oclCoeffAux = (void *)f_coeffAux; 61120cf1dd8SToby Isaac } else { 61220cf1dd8SToby Isaac oclCoeffAux = NULL; 61320cf1dd8SToby Isaac } 61420cf1dd8SToby Isaac oclInvJ = (void *)f_invJ; 61520cf1dd8SToby Isaac oclDetJ = (void *)f_detJ; 6169371c9d4SSatish Balay } break; 6179371c9d4SSatish Balay case PETSC_DOUBLE: { 61820cf1dd8SToby Isaac PetscInt c, b, d; 61920cf1dd8SToby Isaac 6209566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(Ne * N_bt, &d_coeff, Ne, &d_coeffAux, Ne * dim * dim, &d_invJ, Ne, &d_detJ)); 62120cf1dd8SToby Isaac for (c = 0; c < Ne; ++c) { 6227be5e748SToby Isaac d_detJ[c] = (double)cgeom->detJ[c]; 623ad540459SPierre Jolivet for (d = 0; d < dim * dim; ++d) d_invJ[c * dim * dim + d] = (double)cgeom->invJ[c * dim * dim + d]; 624ad540459SPierre Jolivet for (b = 0; b < N_bt; ++b) d_coeff[c * N_bt + b] = (double)coefficients[c * N_bt + b]; 62520cf1dd8SToby Isaac } 62620cf1dd8SToby Isaac if (coefficientsAux) { /* Assume P0 */ 627ad540459SPierre Jolivet for (c = 0; c < Ne; ++c) d_coeffAux[c] = (double)coefficientsAux[c]; 62820cf1dd8SToby Isaac } 62920cf1dd8SToby Isaac oclCoeff = (void *)d_coeff; 63020cf1dd8SToby Isaac if (coefficientsAux) { 63120cf1dd8SToby Isaac oclCoeffAux = (void *)d_coeffAux; 63220cf1dd8SToby Isaac } else { 63320cf1dd8SToby Isaac oclCoeffAux = NULL; 63420cf1dd8SToby Isaac } 63520cf1dd8SToby Isaac oclInvJ = (void *)d_invJ; 63620cf1dd8SToby Isaac oclDetJ = (void *)d_detJ; 6379371c9d4SSatish Balay } break; 638d71ae5a4SJacob Faibussowitsch default: 639d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported PETSc type %d", ocl->realType); 64020cf1dd8SToby Isaac } 64120cf1dd8SToby Isaac } else { 64220cf1dd8SToby Isaac PetscInt c, d; 64320cf1dd8SToby Isaac 6449566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Ne * dim * dim, &r_invJ, Ne, &r_detJ)); 64520cf1dd8SToby Isaac for (c = 0; c < Ne; ++c) { 6467be5e748SToby Isaac r_detJ[c] = cgeom->detJ[c]; 647ad540459SPierre Jolivet for (d = 0; d < dim * dim; ++d) r_invJ[c * dim * dim + d] = cgeom->invJ[c * dim * dim + d]; 64820cf1dd8SToby Isaac } 64920cf1dd8SToby Isaac oclCoeff = (void *)coefficients; 65020cf1dd8SToby Isaac oclCoeffAux = (void *)coefficientsAux; 65120cf1dd8SToby Isaac oclInvJ = (void *)r_invJ; 65220cf1dd8SToby Isaac oclDetJ = (void *)r_detJ; 65320cf1dd8SToby Isaac } 654d0609cedSBarry Smith o_coefficients = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne * N_bt * realSize, oclCoeff, &err); 65520cf1dd8SToby Isaac if (coefficientsAux) { 656d0609cedSBarry Smith o_coefficientsAux = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne * realSize, oclCoeffAux, &err); 65720cf1dd8SToby Isaac } else { 658d0609cedSBarry Smith o_coefficientsAux = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY, Ne * realSize, oclCoeffAux, &err); 65920cf1dd8SToby Isaac } 660d0609cedSBarry Smith o_jacobianInverses = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne * dim * dim * realSize, oclInvJ, &err); 661d0609cedSBarry Smith o_jacobianDeterminants = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne * realSize, oclDetJ, &err); 662d0609cedSBarry Smith o_elemVec = clCreateBuffer(ocl->ctx_id, CL_MEM_WRITE_ONLY, Ne * N_bt * realSize, NULL, &err); 66320cf1dd8SToby Isaac /* Kernel launch */ 6649566063dSJacob Faibussowitsch PetscCall(clSetKernelArg(ocl_kernel, 0, sizeof(cl_int), (void *)&N_cb)); 6659566063dSJacob Faibussowitsch PetscCall(clSetKernelArg(ocl_kernel, 1, sizeof(cl_mem), (void *)&o_coefficients)); 6669566063dSJacob Faibussowitsch PetscCall(clSetKernelArg(ocl_kernel, 2, sizeof(cl_mem), (void *)&o_coefficientsAux)); 6679566063dSJacob Faibussowitsch PetscCall(clSetKernelArg(ocl_kernel, 3, sizeof(cl_mem), (void *)&o_jacobianInverses)); 6689566063dSJacob Faibussowitsch PetscCall(clSetKernelArg(ocl_kernel, 4, sizeof(cl_mem), (void *)&o_jacobianDeterminants)); 6699566063dSJacob Faibussowitsch PetscCall(clSetKernelArg(ocl_kernel, 5, sizeof(cl_mem), (void *)&o_elemVec)); 6709566063dSJacob Faibussowitsch PetscCall(clEnqueueNDRangeKernel(ocl->queue_id, ocl_kernel, 3, NULL, global_work_size, local_work_size, 0, NULL, &ocl_ev)); 67120cf1dd8SToby Isaac /* Read data back from device */ 67220cf1dd8SToby Isaac if (sizeof(PetscReal) != realSize) { 67320cf1dd8SToby Isaac switch (ocl->realType) { 6749371c9d4SSatish Balay case PETSC_FLOAT: { 67520cf1dd8SToby Isaac float *elem; 67620cf1dd8SToby Isaac PetscInt c, b; 67720cf1dd8SToby Isaac 6789566063dSJacob Faibussowitsch PetscCall(PetscFree4(f_coeff, f_coeffAux, f_invJ, f_detJ)); 6799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Ne * N_bt, &elem)); 6809566063dSJacob Faibussowitsch PetscCall(clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne * N_bt * realSize, elem, 0, NULL, NULL)); 68120cf1dd8SToby Isaac for (c = 0; c < Ne; ++c) { 682ad540459SPierre Jolivet for (b = 0; b < N_bt; ++b) elemVec[c * N_bt + b] = (PetscScalar)elem[c * N_bt + b]; 68320cf1dd8SToby Isaac } 6849566063dSJacob Faibussowitsch PetscCall(PetscFree(elem)); 6859371c9d4SSatish Balay } break; 6869371c9d4SSatish Balay case PETSC_DOUBLE: { 68720cf1dd8SToby Isaac double *elem; 68820cf1dd8SToby Isaac PetscInt c, b; 68920cf1dd8SToby Isaac 6909566063dSJacob Faibussowitsch PetscCall(PetscFree4(d_coeff, d_coeffAux, d_invJ, d_detJ)); 6919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Ne * N_bt, &elem)); 6929566063dSJacob Faibussowitsch PetscCall(clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne * N_bt * realSize, elem, 0, NULL, NULL)); 69320cf1dd8SToby Isaac for (c = 0; c < Ne; ++c) { 694ad540459SPierre Jolivet for (b = 0; b < N_bt; ++b) elemVec[c * N_bt + b] = (PetscScalar)elem[c * N_bt + b]; 69520cf1dd8SToby Isaac } 6969566063dSJacob Faibussowitsch PetscCall(PetscFree(elem)); 6979371c9d4SSatish Balay } break; 698d71ae5a4SJacob Faibussowitsch default: 699d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported PETSc type %d", ocl->realType); 70020cf1dd8SToby Isaac } 70120cf1dd8SToby Isaac } else { 7029566063dSJacob Faibussowitsch PetscCall(PetscFree2(r_invJ, r_detJ)); 7039566063dSJacob Faibussowitsch PetscCall(clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne * N_bt * realSize, elemVec, 0, NULL, NULL)); 70420cf1dd8SToby Isaac } 70520cf1dd8SToby Isaac /* Log performance */ 7069566063dSJacob Faibussowitsch PetscCall(clGetEventProfilingInfo(ocl_ev, CL_PROFILING_COMMAND_START, sizeof(cl_ulong), &ns_start, NULL)); 7079566063dSJacob Faibussowitsch PetscCall(clGetEventProfilingInfo(ocl_ev, CL_PROFILING_COMMAND_END, sizeof(cl_ulong), &ns_end, NULL)); 70820cf1dd8SToby Isaac f0Flops = 0; 70920cf1dd8SToby Isaac switch (ocl->op) { 710d71ae5a4SJacob Faibussowitsch case LAPLACIAN: 711d71ae5a4SJacob Faibussowitsch f1Flops = useAux ? dim : 0; 712d71ae5a4SJacob Faibussowitsch break; 713d71ae5a4SJacob Faibussowitsch case ELASTICITY: 714d71ae5a4SJacob Faibussowitsch f1Flops = 2 * dim * dim; 715d71ae5a4SJacob Faibussowitsch break; 71620cf1dd8SToby Isaac } 7179371c9d4SSatish Balay numFlops = Ne * (N_q * (N_b * N_comp * ((useField ? 2 : 0) + (useFieldDer ? 2 * dim * (dim + 1) : 0)) 71820cf1dd8SToby Isaac /*+ 71920cf1dd8SToby Isaac N_ba*N_compa*((useFieldAux ? 2 : 0) + (useFieldDerAux ? 2*dim*(dim + 1) : 0))*/ 7209371c9d4SSatish Balay + N_comp * ((useF0 ? f0Flops + 2 : 0) + (useF1 ? f1Flops + 2 * dim : 0))) + 72120cf1dd8SToby Isaac N_b * ((useF0 ? 2 : 0) + (useF1 ? 2 * dim * (dim + 1) : 0))); 7229566063dSJacob Faibussowitsch PetscCall(PetscFEOpenCLLogResidual(fem, (ns_end - ns_start) * 1.0e-9, numFlops)); 72320cf1dd8SToby Isaac /* Cleanup */ 7249566063dSJacob Faibussowitsch PetscCall(clReleaseMemObject(o_coefficients)); 7259566063dSJacob Faibussowitsch PetscCall(clReleaseMemObject(o_coefficientsAux)); 7269566063dSJacob Faibussowitsch PetscCall(clReleaseMemObject(o_jacobianInverses)); 7279566063dSJacob Faibussowitsch PetscCall(clReleaseMemObject(o_jacobianDeterminants)); 7289566063dSJacob Faibussowitsch PetscCall(clReleaseMemObject(o_elemVec)); 7299566063dSJacob Faibussowitsch PetscCall(clReleaseKernel(ocl_kernel)); 7309566063dSJacob Faibussowitsch PetscCall(clReleaseProgram(ocl_prog)); 73120cf1dd8SToby Isaac PetscFunctionReturn(0); 73220cf1dd8SToby Isaac } 73320cf1dd8SToby Isaac 734526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE); 735526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation); 7362e47ffbbSToby Isaac 737d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_OpenCL(PetscFE fem) 738d71ae5a4SJacob Faibussowitsch { 73920cf1dd8SToby Isaac PetscFunctionBegin; 74020cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 74120cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 74220cf1dd8SToby Isaac fem->ops->view = NULL; 74320cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_OpenCL; 74420cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 745ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 74620cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_OpenCL; 74720cf1dd8SToby Isaac fem->ops->integratebdresidual = NULL /* PetscFEIntegrateBdResidual_OpenCL */; 74820cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_OpenCL */; 74920cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 75020cf1dd8SToby Isaac PetscFunctionReturn(0); 75120cf1dd8SToby Isaac } 75220cf1dd8SToby Isaac 75320cf1dd8SToby Isaac /*MC 754*dce8aebaSBarry Smith PETSCFEOPENCL = "opencl" - A `PetscFEType` that integrates using a vectorized OpenCL implementation 75520cf1dd8SToby Isaac 75620cf1dd8SToby Isaac Level: intermediate 75720cf1dd8SToby Isaac 758db781477SPatrick Sanan .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 75920cf1dd8SToby Isaac M*/ 76020cf1dd8SToby Isaac 761d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_OpenCL(PetscFE fem) 762d71ae5a4SJacob Faibussowitsch { 76320cf1dd8SToby Isaac PetscFE_OpenCL *ocl; 76420cf1dd8SToby Isaac cl_uint num_platforms; 76520cf1dd8SToby Isaac cl_platform_id platform_ids[42]; 76620cf1dd8SToby Isaac cl_uint num_devices; 76720cf1dd8SToby Isaac cl_device_id device_ids[42]; 7682da392ccSBarry Smith cl_int err; 76920cf1dd8SToby Isaac 77020cf1dd8SToby Isaac PetscFunctionBegin; 77120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 7724dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ocl)); 77320cf1dd8SToby Isaac fem->data = ocl; 77420cf1dd8SToby Isaac 77520cf1dd8SToby Isaac /* Init Platform */ 7769566063dSJacob Faibussowitsch PetscCall(clGetPlatformIDs(42, platform_ids, &num_platforms)); 77728b400f6SJacob Faibussowitsch PetscCheck(num_platforms, PetscObjectComm((PetscObject)fem), PETSC_ERR_SUP, "No OpenCL platform found."); 77820cf1dd8SToby Isaac ocl->pf_id = platform_ids[0]; 77920cf1dd8SToby Isaac /* Init Device */ 7809566063dSJacob Faibussowitsch PetscCall(clGetDeviceIDs(ocl->pf_id, CL_DEVICE_TYPE_ALL, 42, device_ids, &num_devices)); 78128b400f6SJacob Faibussowitsch PetscCheck(num_devices, PetscObjectComm((PetscObject)fem), PETSC_ERR_SUP, "No OpenCL device found."); 78220cf1dd8SToby Isaac ocl->dev_id = device_ids[0]; 78320cf1dd8SToby Isaac /* Create context with one command queue */ 7849371c9d4SSatish Balay ocl->ctx_id = clCreateContext(0, 1, &(ocl->dev_id), NULL, NULL, &err); 7859371c9d4SSatish Balay PetscCall(err); 7869371c9d4SSatish Balay ocl->queue_id = clCreateCommandQueue(ocl->ctx_id, ocl->dev_id, CL_QUEUE_PROFILING_ENABLE, &err); 7879371c9d4SSatish Balay PetscCall(err); 78820cf1dd8SToby Isaac /* Types */ 78920cf1dd8SToby Isaac ocl->realType = PETSC_FLOAT; 79020cf1dd8SToby Isaac /* Register events */ 7919566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("OpenCL FEResidual", PETSCFE_CLASSID, &ocl->residualEvent)); 79220cf1dd8SToby Isaac /* Equation handling */ 79320cf1dd8SToby Isaac ocl->op = LAPLACIAN; 79420cf1dd8SToby Isaac 7959566063dSJacob Faibussowitsch PetscCall(PetscFEInitialize_OpenCL(fem)); 79620cf1dd8SToby Isaac PetscFunctionReturn(0); 79720cf1dd8SToby Isaac } 79820cf1dd8SToby Isaac 7992b99622eSMatthew G. Knepley /*@ 800*dce8aebaSBarry Smith PetscFEOpenCLSetRealType - Set the scalar type for running on the OpenCL accelerator 8012b99622eSMatthew G. Knepley 8022b99622eSMatthew G. Knepley Input Parameters: 803*dce8aebaSBarry Smith + fem - The `PetscFE` 8042b99622eSMatthew G. Knepley - realType - The scalar type 8052b99622eSMatthew G. Knepley 8062b99622eSMatthew G. Knepley Level: developer 8072b99622eSMatthew G. Knepley 808*dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEOpenCLGetRealType()` 8092b99622eSMatthew G. Knepley @*/ 810d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEOpenCLSetRealType(PetscFE fem, PetscDataType realType) 811d71ae5a4SJacob Faibussowitsch { 81220cf1dd8SToby Isaac PetscFE_OpenCL *ocl = (PetscFE_OpenCL *)fem->data; 81320cf1dd8SToby Isaac 81420cf1dd8SToby Isaac PetscFunctionBegin; 81520cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 81620cf1dd8SToby Isaac ocl->realType = realType; 81720cf1dd8SToby Isaac PetscFunctionReturn(0); 81820cf1dd8SToby Isaac } 81920cf1dd8SToby Isaac 8202b99622eSMatthew G. Knepley /*@ 821*dce8aebaSBarry Smith PetscFEOpenCLGetRealType - Get the scalar type for running on the OpenCL accelerator 8222b99622eSMatthew G. Knepley 8232b99622eSMatthew G. Knepley Input Parameter: 824*dce8aebaSBarry Smith . fem - The `PetscFE` 8252b99622eSMatthew G. Knepley 8262b99622eSMatthew G. Knepley Output Parameter: 8272b99622eSMatthew G. Knepley . realType - The scalar type 8282b99622eSMatthew G. Knepley 8292b99622eSMatthew G. Knepley Level: developer 8302b99622eSMatthew G. Knepley 831*dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEOpenCLSetRealType()` 8322b99622eSMatthew G. Knepley @*/ 833d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEOpenCLGetRealType(PetscFE fem, PetscDataType *realType) 834d71ae5a4SJacob Faibussowitsch { 83520cf1dd8SToby Isaac PetscFE_OpenCL *ocl = (PetscFE_OpenCL *)fem->data; 83620cf1dd8SToby Isaac 83720cf1dd8SToby Isaac PetscFunctionBegin; 83820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 83920cf1dd8SToby Isaac PetscValidPointer(realType, 2); 84020cf1dd8SToby Isaac *realType = ocl->realType; 84120cf1dd8SToby Isaac PetscFunctionReturn(0); 84220cf1dd8SToby Isaac } 84320cf1dd8SToby Isaac 84420cf1dd8SToby Isaac #endif /* PETSC_HAVE_OPENCL */ 845