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