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 /// @file 18 /// Setup DM for Navier-Stokes example using PETSc 19 20 #include "../navierstokes.h" 21 22 // Create mesh 23 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, DM *dm) { 24 PetscErrorCode ierr; 25 PetscFunctionBeginUser; 26 // Create DMPLEX 27 ierr = DMCreate(comm, dm); CHKERRQ(ierr); 28 ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr); 29 // Set Tensor elements 30 ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr); 31 // Set CL options 32 ierr = DMSetFromOptions(*dm); CHKERRQ(ierr); 33 ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 34 PetscFunctionReturn(0); 35 } 36 37 // Setup DM 38 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, 39 SimpleBC bc, Physics phys, void *setup_ctx) { 40 PetscErrorCode ierr; 41 PetscFunctionBeginUser; 42 { 43 // Configure the finite element space and boundary conditions 44 PetscFE fe; 45 PetscInt num_comp_q = 5; 46 ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, 47 PETSC_FALSE, degree, PETSC_DECIDE, 48 &fe); CHKERRQ(ierr); 49 ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 50 ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 51 ierr = DMCreateDS(dm); CHKERRQ(ierr); 52 { 53 /* create FE field for coordinates */ 54 PetscFE fe_coords; 55 PetscInt num_comp_coord; 56 ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 57 ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_coord, 58 PETSC_FALSE, 1, 1, &fe_coords); CHKERRQ(ierr); 59 ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr); 60 ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr); 61 } 62 ierr = problem->bc_func(dm, bc, phys, setup_ctx); 63 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 64 CHKERRQ(ierr); 65 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 66 } 67 { 68 // Empty name for conserved field (because there is only one field) 69 PetscSection section; 70 ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 71 ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 72 ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 73 CHKERRQ(ierr); 74 ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X"); 75 CHKERRQ(ierr); 76 ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y"); 77 CHKERRQ(ierr); 78 ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z"); 79 CHKERRQ(ierr); 80 ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density"); 81 CHKERRQ(ierr); 82 } 83 PetscFunctionReturn(0); 84 } 85 86 // Refine DM for high-order viz 87 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, 88 SimpleBC bc, Physics phys, void *setup_ctx) { 89 PetscErrorCode ierr; 90 DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 91 VecType vec_type; 92 PetscFunctionBeginUser; 93 94 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 95 96 dm_hierarchy[0] = dm; 97 for (PetscInt i = 0, d = user->app_ctx->degree; 98 i < user->app_ctx->viz_refine; i++) { 99 Mat interp_next; 100 ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]); 101 CHKERRQ(ierr); 102 ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr); 103 ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr); 104 ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr); 105 d = (d + 1) / 2; 106 if (i + 1 == user->app_ctx->viz_refine) d = 1; 107 ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 108 ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr); 109 ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx); 110 CHKERRQ(ierr); 111 ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next, 112 NULL); CHKERRQ(ierr); 113 if (!i) user->interp_viz = interp_next; 114 else { 115 Mat C; 116 ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, 117 PETSC_DECIDE, &C); CHKERRQ(ierr); 118 ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 119 ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr); 120 user->interp_viz = C; 121 } 122 } 123 for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) { 124 ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr); 125 } 126 user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 127 128 PetscFunctionReturn(0); 129 } 130