177841947SLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 277841947SLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 377841947SLeila Ghaffari // reserved. See files LICENSE and NOTICE for details. 477841947SLeila Ghaffari // 577841947SLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software 677841947SLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral 777841947SLeila Ghaffari // element discretizations for exascale applications. For more information and 877841947SLeila Ghaffari // source code availability see http://github.com/ceed. 977841947SLeila Ghaffari // 1077841947SLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 1177841947SLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office 1277841947SLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for 1377841947SLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including 1477841947SLeila Ghaffari // software, applications, hardware, advanced system engineering and early 1577841947SLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative. 1677841947SLeila Ghaffari 1777841947SLeila Ghaffari /// @file 1877841947SLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc 1977841947SLeila Ghaffari 2077841947SLeila Ghaffari #include "../navierstokes.h" 2177841947SLeila Ghaffari 2277841947SLeila Ghaffari // Read mesh and distribute DM in parallel 2377841947SLeila Ghaffari PetscErrorCode CreateDistributedDM(MPI_Comm comm, ProblemData *problem, 2477841947SLeila Ghaffari SetupContext setup_ctx, DM *dm) { 2577841947SLeila Ghaffari DM dist_mesh = NULL; 2677841947SLeila Ghaffari PetscPartitioner part; 27b868981dSJed Brown PetscInt dim = problem->dim, faces[3] = {3, 3, 3}; 2877841947SLeila Ghaffari const PetscReal scale[3] = {setup_ctx->lx, setup_ctx->ly, setup_ctx->lz}; 2977841947SLeila Ghaffari PetscErrorCode ierr; 3077841947SLeila Ghaffari PetscFunctionBeginUser; 3177841947SLeila Ghaffari 32b868981dSJed Brown ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", 33b868981dSJed Brown faces, &dim, NULL); CHKERRQ(ierr); 344a2566fdSLeila Ghaffari if (!dim) dim = problem->dim; 35b868981dSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, faces, NULL, scale, 3677841947SLeila Ghaffari NULL, PETSC_TRUE, dm); CHKERRQ(ierr); 3777841947SLeila Ghaffari 3877841947SLeila Ghaffari // Distribute DM in parallel 3977841947SLeila Ghaffari ierr = DMPlexGetPartitioner(*dm, &part); CHKERRQ(ierr); 4077841947SLeila Ghaffari ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 4177841947SLeila Ghaffari ierr = DMPlexDistribute(*dm, 0, NULL, &dist_mesh); CHKERRQ(ierr); 4277841947SLeila Ghaffari if (dist_mesh) { 4377841947SLeila Ghaffari ierr = DMDestroy(dm); CHKERRQ(ierr); 4477841947SLeila Ghaffari *dm = dist_mesh; 4577841947SLeila Ghaffari } 4677841947SLeila Ghaffari ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 4777841947SLeila Ghaffari 4877841947SLeila Ghaffari PetscFunctionReturn(0); 4977841947SLeila Ghaffari } 5077841947SLeila Ghaffari 5177841947SLeila Ghaffari // Setup DM 5277841947SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, 5377841947SLeila Ghaffari SimpleBC bc, Physics phys, void *setup_ctx) { 5477841947SLeila Ghaffari PetscErrorCode ierr; 5577841947SLeila Ghaffari PetscFunctionBeginUser; 5677841947SLeila Ghaffari { 5777841947SLeila Ghaffari // Configure the finite element space and boundary conditions 5877841947SLeila Ghaffari PetscFE fe; 5977841947SLeila Ghaffari PetscInt num_comp_q = 5; 6077841947SLeila Ghaffari ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, 6177841947SLeila Ghaffari PETSC_FALSE, degree, PETSC_DECIDE, 6277841947SLeila Ghaffari &fe); CHKERRQ(ierr); 6377841947SLeila Ghaffari ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 6477841947SLeila Ghaffari ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 6577841947SLeila Ghaffari ierr = DMCreateDS(dm); CHKERRQ(ierr); 66*7ed3e4cdSJeremy L Thompson { 67*7ed3e4cdSJeremy L Thompson /* create FE field for coordinates */ 68*7ed3e4cdSJeremy L Thompson PetscFE fe_coords; 69*7ed3e4cdSJeremy L Thompson PetscInt num_comp_coord; 70*7ed3e4cdSJeremy L Thompson ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 71*7ed3e4cdSJeremy L Thompson ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_coord, 72*7ed3e4cdSJeremy L Thompson PETSC_FALSE, 1, 1, &fe_coords); CHKERRQ(ierr); 73*7ed3e4cdSJeremy L Thompson ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr); 74*7ed3e4cdSJeremy L Thompson ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr); 75*7ed3e4cdSJeremy L Thompson } 7677841947SLeila Ghaffari ierr = problem->bc_func(dm, bc, phys, setup_ctx); 7777841947SLeila Ghaffari ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 7877841947SLeila Ghaffari CHKERRQ(ierr); 7977841947SLeila Ghaffari ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 8077841947SLeila Ghaffari } 8177841947SLeila Ghaffari { 8277841947SLeila Ghaffari // Empty name for conserved field (because there is only one field) 8377841947SLeila Ghaffari PetscSection section; 8477841947SLeila Ghaffari ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 8577841947SLeila Ghaffari ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 8677841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 8777841947SLeila Ghaffari CHKERRQ(ierr); 8877841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX"); 8977841947SLeila Ghaffari CHKERRQ(ierr); 9077841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY"); 9177841947SLeila Ghaffari CHKERRQ(ierr); 9277841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ"); 9377841947SLeila Ghaffari CHKERRQ(ierr); 9477841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity"); 9577841947SLeila Ghaffari CHKERRQ(ierr); 9677841947SLeila Ghaffari } 9777841947SLeila Ghaffari PetscFunctionReturn(0); 9877841947SLeila Ghaffari } 9977841947SLeila Ghaffari 10077841947SLeila Ghaffari // Refine DM for high-order viz 10177841947SLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, 10277841947SLeila Ghaffari SimpleBC bc, Physics phys, void *setup_ctx) { 10377841947SLeila Ghaffari PetscErrorCode ierr; 10477841947SLeila Ghaffari DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 10577841947SLeila Ghaffari VecType vec_type; 10677841947SLeila Ghaffari PetscFunctionBeginUser; 10777841947SLeila Ghaffari 10877841947SLeila Ghaffari ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 10977841947SLeila Ghaffari 11077841947SLeila Ghaffari dm_hierarchy[0] = dm; 11177841947SLeila Ghaffari for (PetscInt i = 0, d = user->app_ctx->degree; 11277841947SLeila Ghaffari i < user->app_ctx->viz_refine; i++) { 11377841947SLeila Ghaffari Mat interp_next; 11477841947SLeila Ghaffari ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]); 11577841947SLeila Ghaffari CHKERRQ(ierr); 11677841947SLeila Ghaffari ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr); 11777841947SLeila Ghaffari ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr); 11877841947SLeila Ghaffari ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr); 11977841947SLeila Ghaffari d = (d + 1) / 2; 12077841947SLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 12177841947SLeila Ghaffari ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 12277841947SLeila Ghaffari ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr); 12377841947SLeila Ghaffari ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx); 12477841947SLeila Ghaffari CHKERRQ(ierr); 12577841947SLeila Ghaffari ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next, 12677841947SLeila Ghaffari NULL); CHKERRQ(ierr); 12777841947SLeila Ghaffari if (!i) user->interp_viz = interp_next; 12877841947SLeila Ghaffari else { 12977841947SLeila Ghaffari Mat C; 13077841947SLeila Ghaffari ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, 13177841947SLeila Ghaffari PETSC_DECIDE, &C); CHKERRQ(ierr); 13277841947SLeila Ghaffari ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 13377841947SLeila Ghaffari ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr); 13477841947SLeila Ghaffari user->interp_viz = C; 13577841947SLeila Ghaffari } 13677841947SLeila Ghaffari } 13777841947SLeila Ghaffari for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) { 13877841947SLeila Ghaffari ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr); 13977841947SLeila Ghaffari } 14077841947SLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 14177841947SLeila Ghaffari 14277841947SLeila Ghaffari PetscFunctionReturn(0); 14377841947SLeila Ghaffari } 144