1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 // libCEED + PETSc Example: Surface Area 9 // 10 // This example demonstrates a simple usage of libCEED with PETSc to calculate the surface area of a simple closed surface, such as the one of a cube 11 // or a tensor-product discrete sphere via the mass operator. 12 // 13 // The code uses higher level communication protocols in DMPlex. 14 // 15 // Build with: 16 // 17 // make area [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 18 // 19 // Sample runs: 20 // Sequential: 21 // 22 // ./area -problem cube -degree 3 -dm_refine 2 23 // ./area -problem sphere -degree 3 -dm_refine 2 24 // 25 // In parallel: 26 // 27 // mpiexec -n 4 ./area -problem cube -degree 3 -dm_refine 2 28 // mpiexec -n 4 ./area -problem sphere -degree 3 -dm_refine 2 29 // 30 // The above example runs use 2 levels of refinement for the mesh. 31 // Use -dm_refine k, for k levels of uniform refinement. 32 // 33 //TESTARGS -ceed {ceed_resource} -test -degree 3 -dm_refine 1 34 35 /// @file 36 /// libCEED example using the mass operator to compute a cube or a cubed-sphere surface area using PETSc with DMPlex 37 static const char help[] = "Compute surface area of a cube or a cubed-sphere using DMPlex in PETSc\n"; 38 39 #include "area.h" 40 41 #include <ceed.h> 42 #include <petscdmplex.h> 43 #include <petscksp.h> 44 #include <stdbool.h> 45 #include <string.h> 46 47 #include "include/areaproblemdata.h" 48 #include "include/libceedsetup.h" 49 #include "include/matops.h" 50 #include "include/petscutils.h" 51 #include "include/petscversion.h" 52 #include "include/structs.h" 53 54 #ifndef M_PI 55 #define M_PI 3.14159265358979323846 56 #endif 57 58 int main(int argc, char **argv) { 59 MPI_Comm comm; 60 char filename[PETSC_MAX_PATH_LEN], ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 61 PetscInt l_size, g_size, xl_size, 62 q_extra = 1, // default number of extra quadrature points 63 num_comp_x = 3, // number of components of 3D physical coordinates 64 num_comp_u = 1, // dimension of field to which apply mass operator 65 topo_dim = 2, // topological dimension of manifold 66 degree = 3; // default degree for finite element bases 67 PetscBool read_mesh = PETSC_FALSE, test_mode = PETSC_FALSE, simplex = PETSC_FALSE; 68 Vec U, U_loc, V, V_loc; 69 DM dm; 70 OperatorApplyContext op_apply_ctx; 71 Ceed ceed; 72 CeedData ceed_data; 73 ProblemType problem_choice; 74 VecType vec_type; 75 PetscMemType mem_type; 76 77 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 78 comm = PETSC_COMM_WORLD; 79 80 // Read command line options 81 PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc", NULL); 82 problem_choice = SPHERE; 83 PetscCall(PetscOptionsEnum("-problem", "Problem to solve", NULL, problem_types, (PetscEnum)problem_choice, (PetscEnum *)&problem_choice, NULL)); 84 PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL)); 85 PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL)); 86 PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL)); 87 PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh)); 88 PetscCall(PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", NULL, simplex, &simplex, NULL)); 89 PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, °ree, NULL)); 90 PetscOptionsEnd(); 91 92 // Setup DM 93 if (read_mesh) { 94 PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm)); 95 } else { 96 // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box 97 PetscCall(DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm)); 98 if (problem_choice == CUBE) { 99 PetscCall(DMPlexCreateCoordinateSpace(dm, 1, PETSC_TRUE, NULL)); 100 } 101 // Set the object name 102 PetscCall(PetscObjectSetName((PetscObject)dm, problem_types[problem_choice])); 103 // Refine DMPlex with uniform refinement using runtime option -dm_refine 104 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 105 PetscCall(DMSetFromOptions(dm)); 106 // View DMPlex via runtime option 107 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 108 } 109 110 // Create DM 111 PetscCall(SetupDMByDegree(dm, degree, q_extra, num_comp_u, topo_dim, false)); 112 113 // Create vectors 114 PetscCall(DMCreateGlobalVector(dm, &U)); 115 PetscCall(VecGetLocalSize(U, &l_size)); 116 PetscCall(VecGetSize(U, &g_size)); 117 PetscCall(DMCreateLocalVector(dm, &U_loc)); 118 PetscCall(VecGetSize(U_loc, &xl_size)); 119 PetscCall(VecDuplicate(U, &V)); 120 PetscCall(VecDuplicate(U_loc, &V_loc)); 121 122 // Setup op_apply_ctx structure 123 PetscCall(PetscMalloc1(1, &op_apply_ctx)); 124 125 // Set up libCEED 126 CeedInit(ceed_resource, &ceed); 127 CeedMemType mem_type_backend; 128 CeedGetPreferredMemType(ceed, &mem_type_backend); 129 130 PetscCall(DMGetVecType(dm, &vec_type)); 131 if (!vec_type) { // Not yet set by op_apply_ctx -dm_vec_type 132 switch (mem_type_backend) { 133 case CEED_MEM_HOST: 134 vec_type = VECSTANDARD; 135 break; 136 case CEED_MEM_DEVICE: { 137 const char *resolved; 138 CeedGetResource(ceed, &resolved); 139 if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 140 else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 141 else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 142 else vec_type = VECSTANDARD; 143 } 144 } 145 PetscCall(DMSetVecType(dm, vec_type)); 146 } 147 148 // Print summary 149 if (!test_mode) { 150 PetscInt P = degree + 1, Q = P + q_extra; 151 const char *used_resource; 152 CeedGetResource(ceed, &used_resource); 153 PetscCall(PetscPrintf(comm, 154 "\n-- libCEED + PETSc Surface Area of a Manifold --\n" 155 " libCEED:\n" 156 " libCEED Backend : %s\n" 157 " libCEED Backend MemType : %s\n" 158 " Mesh:\n" 159 " Solution Order (P) : %" PetscInt_FMT "\n" 160 " Quadrature Order (Q) : %" PetscInt_FMT "\n" 161 " Additional quadrature points (q_extra) : %" PetscInt_FMT "\n" 162 " Global nodes : %" PetscInt_FMT "\n" 163 " DoF per node : %" PetscInt_FMT "\n" 164 " Global DoFs : %" PetscInt_FMT "\n", 165 used_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size / num_comp_u, num_comp_u, g_size)); 166 } 167 168 // Setup libCEED's objects and apply setup operator 169 PetscCall(PetscMalloc1(1, &ceed_data)); 170 PetscCall(SetupLibceedByDegree(dm, ceed, degree, topo_dim, q_extra, num_comp_x, num_comp_u, g_size, xl_size, problem_options[problem_choice], 171 ceed_data, false, (CeedVector)NULL, (CeedVector *)NULL)); 172 173 // Setup output vector 174 PetscCall(VecZeroEntries(V_loc)); 175 PetscCall(VecP2C(V_loc, &mem_type, ceed_data->y_ceed)); 176 177 // Compute the mesh volume using the mass operator: area = 1^T \cdot M \cdot 1 178 if (!test_mode) { 179 PetscCall(PetscPrintf(comm, "Computing the mesh area using the formula: area = 1^T M 1\n")); 180 } 181 182 // Initialize u with ones 183 CeedVectorSetValue(ceed_data->x_ceed, 1.0); 184 185 // Apply the mass operator: 'u' -> 'v' 186 CeedOperatorApply(ceed_data->op_apply, ceed_data->x_ceed, ceed_data->y_ceed, CEED_REQUEST_IMMEDIATE); 187 188 // Gather output vector 189 PetscCall(VecC2P(ceed_data->y_ceed, mem_type, V_loc)); 190 PetscCall(VecZeroEntries(V)); 191 PetscCall(DMLocalToGlobalBegin(dm, V_loc, ADD_VALUES, V)); 192 PetscCall(DMLocalToGlobalEnd(dm, V_loc, ADD_VALUES, V)); 193 194 // Compute and print the sum of the entries of 'v' giving the mesh surface area 195 PetscScalar area; 196 PetscCall(VecSum(V, &area)); 197 198 // Compute the exact surface area and print the result 199 CeedScalar exact_surface_area = 4 * M_PI; 200 if (problem_choice == CUBE) { 201 exact_surface_area = 6 * 2 * 2; // surface of [-1, 1]^3 202 } 203 204 PetscReal error = fabs(area - exact_surface_area); 205 PetscReal tol = 5e-6; 206 if (!test_mode || error > tol) { 207 PetscCall(PetscPrintf(comm, "Exact mesh surface area : % .14g\n", exact_surface_area)); 208 PetscCall(PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area)); 209 PetscCall(PetscPrintf(comm, "Area error : % .14g\n", error)); 210 } 211 212 // Cleanup 213 PetscCall(DMDestroy(&dm)); 214 PetscCall(VecDestroy(&U)); 215 PetscCall(VecDestroy(&U_loc)); 216 PetscCall(VecDestroy(&V)); 217 PetscCall(VecDestroy(&V_loc)); 218 PetscCall(PetscFree(op_apply_ctx)); 219 PetscCall(CeedDataDestroy(0, ceed_data)); 220 CeedDestroy(&ceed); 221 return PetscFinalize(); 222 } 223