1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 // libCEED + PETSc Example: Surface Area 18 // 19 // This example demonstrates a simple usage of libCEED with PETSc to calculate 20 // the surface area of a simple closed surface, such as the one of a cube or a 21 // tensor-product discrete sphere via the mass operator. 22 // 23 // The code uses higher level communication protocols in DMPlex. 24 // 25 // Build with: 26 // 27 // make area [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 28 // 29 // Sample runs: 30 // Sequential: 31 // 32 // ./area -problem cube -degree 3 -dm_refine 2 33 // ./area -problem sphere -degree 3 -dm_refine 2 34 // 35 // In parallel: 36 // 37 // mpiexec -n 4 ./area -problem cube -degree 3 -dm_refine 2 38 // mpiexec -n 4 ./area -problem sphere -degree 3 -dm_refine 2 39 // 40 // The above example runs use 2 levels of refinement for the mesh. 41 // Use -dm_refine k, for k levels of uniform refinement. 42 // 43 //TESTARGS -ceed {ceed_resource} -test -degree 3 -dm_refine 1 44 45 /// @file 46 /// libCEED example using the mass operator to compute a cube or a cubed-sphere surface area using PETSc with DMPlex 47 static const char help[] = 48 "Compute surface area of a cube or a cubed-sphere using DMPlex in PETSc\n"; 49 50 #include <stdbool.h> 51 #include <string.h> 52 #include <ceed.h> 53 #include <petsc.h> 54 #include <petscdmplex.h> 55 56 #include "area.h" 57 #include "include/areaproblemdata.h" 58 #include "include/petscmacros.h" 59 #include "include/petscutils.h" 60 #include "include/matops.h" 61 #include "include/structs.h" 62 #include "include/libceedsetup.h" 63 64 #if PETSC_VERSION_LT(3,12,0) 65 #ifdef PETSC_HAVE_CUDA 66 #include <petsccuda.h> 67 // Note: With PETSc prior to version 3.12.0, providing the source path to 68 // include 'cublas_v2.h' will be needed to use 'petsccuda.h'. 69 #endif 70 #endif 71 72 #ifndef M_PI 73 # define M_PI 3.14159265358979323846 74 #endif 75 76 int main(int argc, char **argv) { 77 PetscInt ierr; 78 MPI_Comm comm; 79 char filename[PETSC_MAX_PATH_LEN], 80 ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 81 PetscInt l_size, g_size, xl_size, 82 q_extra = 1, // default number of extra quadrature points 83 num_comp_x = 3, // number of components of 3D physical coordinates 84 num_comp_u = 1, // dimension of field to which apply mass operator 85 topo_dim = 2, // topological dimension of manifold 86 degree = 3; // default degree for finite element bases 87 PetscBool read_mesh = PETSC_FALSE, 88 test_mode = PETSC_FALSE, 89 simplex = PETSC_FALSE; 90 Vec U, U_loc, V, V_loc; 91 DM dm; 92 UserO user; 93 Ceed ceed; 94 CeedData ceed_data; 95 ProblemType problem_choice; 96 VecType vec_type; 97 PetscMemType mem_type; 98 99 ierr = PetscInitialize(&argc, &argv, NULL, help); 100 if (ierr) return ierr; 101 comm = PETSC_COMM_WORLD; 102 103 // Read command line options 104 ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc", 105 NULL); 106 CHKERRQ(ierr); 107 problem_choice = SPHERE; 108 ierr = PetscOptionsEnum("-problem", 109 "Problem to solve", NULL, 110 problem_types, (PetscEnum)problem_choice, 111 (PetscEnum *)&problem_choice, 112 NULL); CHKERRQ(ierr); 113 ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points", 114 NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr); 115 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 116 NULL, ceed_resource, ceed_resource, 117 sizeof(ceed_resource), NULL); CHKERRQ(ierr); 118 ierr = PetscOptionsBool("-test", 119 "Testing mode (do not print unless error is large)", 120 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 121 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 122 filename, filename, sizeof(filename), &read_mesh); 123 CHKERRQ(ierr); 124 ierr = PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", 125 NULL, simplex, &simplex, NULL); CHKERRQ(ierr); 126 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 127 NULL, degree, °ree, NULL); CHKERRQ(ierr); 128 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 129 130 // Setup DM 131 if (read_mesh) { 132 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm); 133 CHKERRQ(ierr); 134 } else { 135 // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box 136 ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm); 137 CHKERRQ(ierr); 138 // Set the object name 139 ierr = PetscObjectSetName((PetscObject)dm, problem_types[problem_choice]); 140 CHKERRQ(ierr); 141 // Distribute mesh over processes 142 { 143 DM dm_dist = NULL; 144 PetscPartitioner part; 145 146 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 147 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 148 ierr = DMPlexDistribute(dm, 0, NULL, &dm_dist); CHKERRQ(ierr); 149 if (dm_dist) { 150 ierr = DMDestroy(&dm); CHKERRQ(ierr); 151 dm = dm_dist; 152 } 153 } 154 // Refine DMPlex with uniform refinement using runtime option -dm_refine 155 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 156 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 157 if (problem_choice == SPHERE) { 158 ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr); 159 } 160 // View DMPlex via runtime option 161 ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 162 } 163 164 // Create DM 165 ierr = SetupDMByDegree(dm, degree, num_comp_u, topo_dim, false, 166 (BCFunction)NULL); 167 CHKERRQ(ierr); 168 169 // Create vectors 170 ierr = DMCreateGlobalVector(dm, &U); CHKERRQ(ierr); 171 ierr = VecGetLocalSize(U, &l_size); CHKERRQ(ierr); 172 ierr = VecGetSize(U, &g_size); CHKERRQ(ierr); 173 ierr = DMCreateLocalVector(dm, &U_loc); CHKERRQ(ierr); 174 ierr = VecGetSize(U_loc, &xl_size); CHKERRQ(ierr); 175 ierr = VecDuplicate(U, &V); CHKERRQ(ierr); 176 ierr = VecDuplicate(U_loc, &V_loc); CHKERRQ(ierr); 177 178 // Setup user structure 179 ierr = PetscMalloc1(1, &user); CHKERRQ(ierr); 180 181 // Set up libCEED 182 CeedInit(ceed_resource, &ceed); 183 CeedMemType mem_type_backend; 184 CeedGetPreferredMemType(ceed, &mem_type_backend); 185 186 ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 187 if (!vec_type) { // Not yet set by user -dm_vec_type 188 switch (mem_type_backend) { 189 case CEED_MEM_HOST: vec_type = VECSTANDARD; break; 190 case CEED_MEM_DEVICE: { 191 const char *resolved; 192 CeedGetResource(ceed, &resolved); 193 if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 194 else if (strstr(resolved, "/gpu/hip/occa")) 195 vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 196 else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 197 else vec_type = VECSTANDARD; 198 } 199 } 200 ierr = DMSetVecType(dm, vec_type); CHKERRQ(ierr); 201 } 202 203 // Print summary 204 if (!test_mode) { 205 PetscInt P = degree + 1, Q = P + q_extra; 206 const char *used_resource; 207 CeedGetResource(ceed, &used_resource); 208 ierr = PetscPrintf(comm, 209 "\n-- libCEED + PETSc Surface Area of a Manifold --\n" 210 " libCEED:\n" 211 " libCEED Backend : %s\n" 212 " libCEED Backend MemType : %s\n" 213 " Mesh:\n" 214 " Number of 1D Basis Nodes (p) : %d\n" 215 " Number of 1D Quadrature Points (q) : %d\n" 216 " Global nodes : %D\n" 217 " DoF per node : %D\n" 218 " Global DoFs : %D\n", 219 used_resource, CeedMemTypes[mem_type_backend], P, Q, 220 g_size/num_comp_u, num_comp_u, g_size); CHKERRQ(ierr); 221 } 222 223 // Setup libCEED's objects and apply setup operator 224 ierr = PetscMalloc1(1, &ceed_data); CHKERRQ(ierr); 225 ierr = SetupLibceedByDegree(dm, ceed, degree, topo_dim, q_extra, num_comp_x, 226 num_comp_u, 227 g_size, xl_size, problem_options[problem_choice], ceed_data, 228 false, (CeedVector)NULL, (CeedVector *)NULL); 229 CHKERRQ(ierr); 230 231 // Setup output vector 232 PetscScalar *v; 233 ierr = VecZeroEntries(V_loc); CHKERRQ(ierr); 234 ierr = VecGetArrayAndMemType(V_loc, &v, &mem_type); CHKERRQ(ierr); 235 CeedVectorSetArray(ceed_data->y_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, 236 v); 237 238 // Compute the mesh volume using the mass operator: area = 1^T \cdot M \cdot 1 239 if (!test_mode) { 240 ierr = PetscPrintf(comm, 241 "Computing the mesh area using the formula: area = 1^T M 1\n"); 242 CHKERRQ(ierr); 243 } 244 245 // Initialize u with ones 246 CeedVectorSetValue(ceed_data->x_ceed, 1.0); 247 248 // Apply the mass operator: 'u' -> 'v' 249 CeedOperatorApply(ceed_data->op_apply, ceed_data->x_ceed, ceed_data->y_ceed, 250 CEED_REQUEST_IMMEDIATE); 251 252 // Gather output vector 253 CeedVectorTakeArray(ceed_data->y_ceed, CEED_MEM_HOST, NULL); 254 ierr = VecRestoreArrayAndMemType(V_loc, &v); CHKERRQ(ierr); 255 ierr = VecZeroEntries(V); CHKERRQ(ierr); 256 ierr = DMLocalToGlobalBegin(dm, V_loc, ADD_VALUES, V); CHKERRQ(ierr); 257 ierr = DMLocalToGlobalEnd(dm, V_loc, ADD_VALUES, V); CHKERRQ(ierr); 258 259 // Compute and print the sum of the entries of 'v' giving the mesh surface area 260 PetscScalar area; 261 ierr = VecSum(V, &area); CHKERRQ(ierr); 262 263 // Compute the exact surface area and print the result 264 CeedScalar exact_surface_area = 4 * M_PI; 265 if (problem_choice == CUBE) { 266 PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube 267 exact_surface_area = 6 * (2*l) * (2*l); 268 } 269 270 PetscReal error = fabs(area - exact_surface_area); 271 PetscReal tol = 5e-6; 272 if (!test_mode || error > tol) { 273 ierr = PetscPrintf(comm, "Exact mesh surface area : % .14g\n", 274 exact_surface_area); 275 CHKERRQ(ierr); 276 ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area); 277 CHKERRQ(ierr); 278 ierr = PetscPrintf(comm, "Area error : % .14g\n", error); 279 CHKERRQ(ierr); 280 } 281 282 // Cleanup 283 ierr = DMDestroy(&dm); CHKERRQ(ierr); 284 ierr = VecDestroy(&U); CHKERRQ(ierr); 285 ierr = VecDestroy(&U_loc); CHKERRQ(ierr); 286 ierr = VecDestroy(&V); CHKERRQ(ierr); 287 ierr = VecDestroy(&V_loc); CHKERRQ(ierr); 288 ierr = PetscFree(user); CHKERRQ(ierr); 289 ierr = CeedDataDestroy(0, ceed_data); CHKERRQ(ierr); 290 CeedDestroy(&ceed); 291 return PetscFinalize(); 292 } 293