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 // 34 // ./area -problem sphere -degree 3 -dm_refine 2 35 // 36 // In parallel: 37 // 38 // mpiexec -n 4 ./area -problem cube -degree 3 -dm_refine 2 39 // 40 // mpiexec -n 4 ./area -problem sphere -degree 3 -dm_refine 2 41 // 42 // The above example runs use 2 levels of refinement for the mesh. 43 // Use -dm_refine k, for k levels of uniform refinement. 44 // 45 //TESTARGS -ceed {ceed_resource} -test -degree 3 -dm_refine 1 46 47 /// @file 48 /// libCEED example using the mass operator to compute a cube or a cubed-sphere surface area using PETSc with DMPlex 49 static const char help[] = 50 "Compute surface area of a cube or a cubed-sphere using DMPlex in PETSc\n"; 51 52 #include <string.h> 53 #include <petscdmplex.h> 54 #include <ceed.h> 55 #include "setuparea.h" 56 57 #ifndef M_PI 58 # define M_PI 3.14159265358979323846 59 #endif 60 61 int main(int argc, char **argv) { 62 PetscInt ierr; 63 MPI_Comm comm; 64 char filename[PETSC_MAX_PATH_LEN], 65 ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 66 PetscInt lsize, gsize, xlsize, 67 qextra = 1, // default number of extra quadrature points 68 ncompx = 3, // number of components of 3D physical coordinates 69 ncompu = 1, // dimension of field to which apply mass operator 70 topodim = 2, // topological dimension of manifold 71 degree = 3; // default degree for finite element bases 72 PetscBool read_mesh = PETSC_FALSE, 73 test_mode = PETSC_FALSE, 74 simplex = PETSC_FALSE; 75 Vec U, Uloc, V, Vloc; 76 DM dm; 77 User user; 78 Ceed ceed; 79 CeedData ceeddata; 80 problemType problemChoice; 81 82 ierr = PetscInitialize(&argc, &argv, NULL, help); 83 if (ierr) return ierr; 84 comm = PETSC_COMM_WORLD; 85 86 // Read command line options 87 ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc", 88 NULL); 89 CHKERRQ(ierr); 90 problemChoice = SPHERE; 91 ierr = PetscOptionsEnum("-problem", 92 "Problem to solve", NULL, 93 problemTypes, (PetscEnum)problemChoice, 94 (PetscEnum *)&problemChoice, 95 NULL); CHKERRQ(ierr); 96 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 97 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 98 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 99 NULL, ceedresource, ceedresource, 100 sizeof(ceedresource), NULL); CHKERRQ(ierr); 101 ierr = PetscOptionsBool("-test", 102 "Testing mode (do not print unless error is large)", 103 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 104 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 105 filename, filename, sizeof(filename), &read_mesh); 106 CHKERRQ(ierr); 107 ierr = PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", 108 NULL, simplex, &simplex, NULL); CHKERRQ(ierr); 109 ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 110 NULL, degree, °ree, NULL); CHKERRQ(ierr); 111 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 112 113 // Setup DM 114 if (read_mesh) { 115 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm); 116 CHKERRQ(ierr); 117 } else { 118 // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box 119 ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, 1., &dm); 120 CHKERRQ(ierr); 121 // Set the object name 122 ierr = PetscObjectSetName((PetscObject)dm, problemTypes[problemChoice]); 123 CHKERRQ(ierr); 124 // Distribute mesh over processes 125 { 126 DM dmDist = NULL; 127 PetscPartitioner part; 128 129 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 130 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 131 ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 132 if (dmDist) { 133 ierr = DMDestroy(&dm); CHKERRQ(ierr); 134 dm = dmDist; 135 } 136 } 137 // Refine DMPlex with uniform refinement using runtime option -dm_refine 138 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 139 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 140 if (problemChoice == SPHERE) { 141 ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr); 142 } 143 // View DMPlex via runtime option 144 ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 145 } 146 147 // Create DM 148 ierr = SetupDMByDegree(dm, degree, ncompu, topodim); CHKERRQ(ierr); 149 150 // Create vectors 151 ierr = DMCreateGlobalVector(dm, &U); CHKERRQ(ierr); 152 ierr = VecGetLocalSize(U, &lsize); CHKERRQ(ierr); 153 ierr = VecGetSize(U, &gsize); CHKERRQ(ierr); 154 ierr = DMCreateLocalVector(dm, &Uloc); CHKERRQ(ierr); 155 ierr = VecGetSize(Uloc, &xlsize); CHKERRQ(ierr); 156 ierr = VecDuplicate(U, &V); CHKERRQ(ierr); 157 ierr = VecDuplicate(Uloc, &Vloc); CHKERRQ(ierr); 158 159 // Setup user structure 160 ierr = PetscMalloc1(1, &user); CHKERRQ(ierr); 161 162 // Set up libCEED 163 CeedInit(ceedresource, &ceed); 164 165 // Print summary 166 if (!test_mode) { 167 PetscInt P = degree + 1, Q = P + qextra; 168 const char *usedresource; 169 CeedGetResource(ceed, &usedresource); 170 ierr = PetscPrintf(comm, 171 "\n-- libCEED + PETSc Surface Area of a Manifold --\n" 172 " libCEED:\n" 173 " libCEED Backend : %s\n" 174 " Mesh:\n" 175 " Number of 1D Basis Nodes (p) : %d\n" 176 " Number of 1D Quadrature Points (q) : %d\n" 177 " Global nodes : %D\n" 178 " DoF per node : %D\n" 179 " Global DoFs : %D\n", 180 usedresource, P, Q, gsize/ncompu, ncompu, gsize); 181 CHKERRQ(ierr); 182 } 183 184 // Setup libCEED's objects and apply setup operator 185 ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr); 186 ierr = SetupLibceedByDegree(dm, ceed, degree, topodim, qextra, 187 ncompx, ncompu, xlsize, problemChoice, 188 ceeddata); CHKERRQ(ierr); 189 190 // Setup output vector 191 PetscScalar *v; 192 ierr = VecZeroEntries(Vloc); CHKERRQ(ierr); 193 ierr = VecGetArray(Vloc, &v); 194 CeedVectorSetArray(ceeddata->vceed, CEED_MEM_HOST, CEED_USE_POINTER, v); 195 196 // Compute the mesh volume using the mass operator: area = 1^T \cdot M \cdot 1 197 if (!test_mode) { 198 ierr = PetscPrintf(comm, 199 "Computing the mesh area using the formula: area = 1^T M 1\n"); 200 CHKERRQ(ierr); 201 } 202 203 // Initialize u with ones 204 CeedVectorSetValue(ceeddata->uceed, 1.0); 205 206 // Apply the mass operator: 'u' -> 'v' 207 CeedOperatorApply(ceeddata->op_apply, ceeddata->uceed, ceeddata->vceed, 208 CEED_REQUEST_IMMEDIATE); 209 210 // Gather output vector 211 CeedVectorTakeArray(ceeddata->vceed, CEED_MEM_HOST, NULL); 212 ierr = VecRestoreArray(Vloc, &v); CHKERRQ(ierr); 213 ierr = VecZeroEntries(V); CHKERRQ(ierr); 214 ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr); 215 ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr); 216 217 // Compute and print the sum of the entries of 'v' giving the mesh surface area 218 PetscScalar area; 219 ierr = VecSum(V, &area); CHKERRQ(ierr); 220 221 // Compute the exact surface area and print the result 222 CeedScalar exact_surfarea = 4 * M_PI; 223 if (problemChoice == CUBE) { 224 PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube 225 exact_surfarea = 6 * (2*l) * (2*l); 226 } 227 228 PetscReal error = fabs(area - exact_surfarea); 229 PetscReal tol = 5e-6; 230 if (!test_mode || error > tol) { 231 ierr = PetscPrintf(comm, "Exact mesh surface area : % .14g\n", 232 exact_surfarea); CHKERRQ(ierr); 233 ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area); 234 CHKERRQ(ierr); 235 ierr = PetscPrintf(comm, "Area error : % .14g\n", 236 (double)error); CHKERRQ(ierr); 237 } 238 239 // Cleanup 240 ierr = DMDestroy(&dm); CHKERRQ(ierr); 241 ierr = VecDestroy(&U); CHKERRQ(ierr); 242 ierr = VecDestroy(&Uloc); CHKERRQ(ierr); 243 ierr = VecDestroy(&V); CHKERRQ(ierr); 244 ierr = VecDestroy(&Vloc); CHKERRQ(ierr); 245 ierr = CeedDataDestroy(ceeddata); CHKERRQ(ierr); 246 CeedDestroy(&ceed); 247 return PetscFinalize(); 248 } 249