19dbcead6SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 29dbcead6SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 39dbcead6SJeremy L Thompson // 49dbcead6SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 59dbcead6SJeremy L Thompson // 69dbcead6SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 79dbcead6SJeremy L Thompson 89dbcead6SJeremy L Thompson // libCEED + PETSc DMSwarm Example 99dbcead6SJeremy L Thompson // 109dbcead6SJeremy L Thompson // This example demonstrates a simple usage of libCEED with DMSwarm. 119dbcead6SJeremy L Thompson // This example combines elements of PETSc src/impls/dm/swam/tutorials/ex1.c and src/impls/dm/swarm/tests/ex6.c 129dbcead6SJeremy L Thompson // 139dbcead6SJeremy L Thompson // Build with: 149dbcead6SJeremy L Thompson // 159dbcead6SJeremy L Thompson // make dmswarm [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 169dbcead6SJeremy L Thompson // 179dbcead6SJeremy L Thompson // Sample runs: 189dbcead6SJeremy L Thompson // 199dbcead6SJeremy L Thompson // ./dmswarm -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_box_lower -1.0,-1.0,-1.0 -dm_plex_simplex 0 -num_comp 2 -swarm gauss 209dbcead6SJeremy L Thompson // 219dbcead6SJeremy L Thompson //TESTARGS(name="Uniform swarm, CG projection") -ceed {ceed_resource} -test -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_box_lower -1.0,-1.0,-1.0 -dm_plex_simplex 0 -dm_plex_hash_location true -num_comp 2 -swarm uniform -solution_order 3 -points_per_cell 125 229dbcead6SJeremy L Thompson //TESTARGS(name="Gauss swarm, lumped projection") -ceed {ceed_resource} -test -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_box_lower -1.0,-1.0,-1.0 -dm_plex_simplex 0 -dm_plex_hash_location true -num_comp 2 -swarm gauss -ksp_type preonly -pc_type jacobi -pc_jacobi_type rowsum -tolerance 9e-2 239dbcead6SJeremy L Thompson 249dbcead6SJeremy L Thompson /// @file 259dbcead6SJeremy L Thompson /// libCEED example using PETSc with DMSwarm 269dbcead6SJeremy L Thompson const char help[] = "libCEED example using PETSc with DMSwarm\n"; 279dbcead6SJeremy L Thompson 289dbcead6SJeremy L Thompson #include <ceed.h> 299dbcead6SJeremy L Thompson #include <math.h> 309dbcead6SJeremy L Thompson #include <petscdmplex.h> 319dbcead6SJeremy L Thompson #include <petscdmswarm.h> 329dbcead6SJeremy L Thompson #include <petscds.h> 339dbcead6SJeremy L Thompson #include <petscfe.h> 349dbcead6SJeremy L Thompson #include <petscksp.h> 359dbcead6SJeremy L Thompson #include <petsc/private/petscfeimpl.h> /* For interpolation */ 369dbcead6SJeremy L Thompson 379dbcead6SJeremy L Thompson #include "include/petscutils.h" 38d31f425aSJeremy L Thompson #include "include/petscversion.h" 39*b6972d74SZach Atkins #include "include/swarmutils.h" 409dbcead6SJeremy L Thompson 419dbcead6SJeremy L Thompson const char DMSwarmPICField_u[] = "u"; 429dbcead6SJeremy L Thompson 439dbcead6SJeremy L Thompson // Target functions 449dbcead6SJeremy L Thompson typedef enum { TARGET_TANH = 0, TARGET_POLYNOMIAL = 1, TARGET_SPHERE = 2 } TargetType; 459dbcead6SJeremy L Thompson static const char *const target_types[] = {"tanh", "polynomial", "sphere", "TargetType", "TARGET", 0}; 469dbcead6SJeremy L Thompson 479dbcead6SJeremy L Thompson typedef PetscErrorCode (*TargetFunc)(PetscInt dim, const PetscScalar x[]); 489dbcead6SJeremy L Thompson typedef PetscErrorCode (*TargetFuncProj)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx); 499dbcead6SJeremy L Thompson 509dbcead6SJeremy L Thompson PetscScalar EvalU_Tanh(PetscInt dim, const PetscScalar x[]); 519dbcead6SJeremy L Thompson PetscScalar EvalU_Poly(PetscInt dim, const PetscScalar x[]); 529dbcead6SJeremy L Thompson PetscScalar EvalU_Sphere(PetscInt dim, const PetscScalar x[]); 539dbcead6SJeremy L Thompson PetscErrorCode EvalU_Tanh_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx); 549dbcead6SJeremy L Thompson PetscErrorCode EvalU_Poly_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx); 559dbcead6SJeremy L Thompson PetscErrorCode EvalU_Sphere_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx); 569dbcead6SJeremy L Thompson 579dbcead6SJeremy L Thompson // Swarm to mesh and mesh to swarm 589dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Petsc(DM dm_swarm, const char *field, Vec U_mesh); 599dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Ceed(DM dm_swarm, const char *field, Vec U_mesh); 609dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCheckSwarmValues(DM dm_swarm, const char *field, PetscScalar tolerance, TargetFuncProj TrueSolution); 619dbcead6SJeremy L Thompson 629dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 639dbcead6SJeremy L Thompson // main driver 649dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 659dbcead6SJeremy L Thompson int main(int argc, char **argv) { 669dbcead6SJeremy L Thompson MPI_Comm comm; 679dbcead6SJeremy L Thompson char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 689dbcead6SJeremy L Thompson PetscBool test_mode = PETSC_FALSE, view_petsc_swarm = PETSC_FALSE, view_ceed_swarm = PETSC_FALSE; 699dbcead6SJeremy L Thompson PetscInt dim = 3, num_comp = 1, num_points = 1728, num_points_per_cell = 64, mesh_order = 1, solution_order = 2, q_extra = 3; 709dbcead6SJeremy L Thompson PetscScalar tolerance = 1E-3; 719dbcead6SJeremy L Thompson DM dm_mesh, dm_swarm; 729dbcead6SJeremy L Thompson Vec U_mesh; 739dbcead6SJeremy L Thompson DMSwarmCeedContext swarm_ceed_context; 749dbcead6SJeremy L Thompson PointSwarmType point_swarm_type = SWARM_UNIFORM; 759dbcead6SJeremy L Thompson TargetType target_type = TARGET_TANH; 769dbcead6SJeremy L Thompson TargetFuncProj target_function_proj = EvalU_Tanh_proj; 779dbcead6SJeremy L Thompson 789dbcead6SJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 799dbcead6SJeremy L Thompson comm = PETSC_COMM_WORLD; 809dbcead6SJeremy L Thompson 819dbcead6SJeremy L Thompson // Read command line options 829dbcead6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "libCEED example using PETSc with DMSwarm", NULL); 839dbcead6SJeremy L Thompson 849dbcead6SJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL)); 859dbcead6SJeremy L Thompson PetscCall( 869dbcead6SJeremy L Thompson PetscOptionsBool("-u_petsc_swarm_view", "View XDMF of swarm values interpolated by PETSc", NULL, view_petsc_swarm, &view_petsc_swarm, NULL)); 879dbcead6SJeremy L Thompson PetscCall( 889dbcead6SJeremy L Thompson PetscOptionsBool("-u_ceed_swarm_view", "View XDMF of swarm values interpolated by libCEED", NULL, view_ceed_swarm, &view_ceed_swarm, NULL)); 899dbcead6SJeremy L Thompson PetscCall(PetscOptionsEnum("-target", "Target field function", NULL, target_types, (PetscEnum)target_type, (PetscEnum *)&target_type, NULL)); 909dbcead6SJeremy L Thompson PetscCall(PetscOptionsInt("-solution_order", "Order of mesh solution space", NULL, solution_order, &solution_order, NULL)); 919dbcead6SJeremy L Thompson PetscCall(PetscOptionsInt("-mesh_order", "Order of mesh coordinate space", NULL, mesh_order, &mesh_order, NULL)); 929dbcead6SJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL)); 939dbcead6SJeremy L Thompson PetscCall(PetscOptionsInt("-num_comp", "Number of components in solution", NULL, num_comp, &num_comp, NULL)); 949dbcead6SJeremy L Thompson PetscCall(PetscOptionsEnum("-swarm", "Swarm points distribution", NULL, point_swarm_types, (PetscEnum)point_swarm_type, 959dbcead6SJeremy L Thompson (PetscEnum *)&point_swarm_type, NULL)); 969dbcead6SJeremy L Thompson { 979dbcead6SJeremy L Thompson PetscBool user_set_num_points_per_cell = PETSC_FALSE; 989dbcead6SJeremy L Thompson PetscInt dim = 3, num_cells_total = 1; 999dbcead6SJeremy L Thompson PetscInt num_cells[] = {1, 1, 1}; 1009dbcead6SJeremy L Thompson 1019dbcead6SJeremy L Thompson PetscCall(PetscOptionsInt("-points_per_cell", "Total number of swarm points in each cell", NULL, num_points_per_cell, &num_points_per_cell, 1029dbcead6SJeremy L Thompson &user_set_num_points_per_cell)); 1039dbcead6SJeremy L Thompson PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL)); 1049dbcead6SJeremy L Thompson PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &dim, NULL)); 1059dbcead6SJeremy L Thompson 1069dbcead6SJeremy L Thompson num_cells_total = num_cells[0] * num_cells[1] * num_cells[2]; 1079dbcead6SJeremy L Thompson PetscCheck(!user_set_num_points_per_cell || point_swarm_type != SWARM_SINUSOIDAL, comm, PETSC_ERR_USER, 1089dbcead6SJeremy L Thompson "Cannot specify points per cell with sinusoidal points locations"); 1099dbcead6SJeremy L Thompson if (!user_set_num_points_per_cell) { 1109dbcead6SJeremy L Thompson PetscCall(PetscOptionsInt("-points", "Total number of swarm points", NULL, num_points, &num_points, NULL)); 1119dbcead6SJeremy L Thompson num_points_per_cell = PetscCeilInt(num_points, num_cells_total); 1129dbcead6SJeremy L Thompson } 1139dbcead6SJeremy L Thompson if (point_swarm_type != SWARM_SINUSOIDAL) { 1149dbcead6SJeremy L Thompson PetscInt num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0)); 1159dbcead6SJeremy L Thompson 1169dbcead6SJeremy L Thompson num_points_per_cell = 1; 1179dbcead6SJeremy L Thompson for (PetscInt i = 0; i < dim; i++) num_points_per_cell *= num_points_per_cell_1d; 1189dbcead6SJeremy L Thompson } 1199dbcead6SJeremy L Thompson num_points = num_points_per_cell * num_cells_total; 1209dbcead6SJeremy L Thompson } 1219dbcead6SJeremy L Thompson PetscCall(PetscOptionsScalar("-tolerance", "Tolerance for swarm point values and projection relative L2 error", NULL, tolerance, &tolerance, NULL)); 1229dbcead6SJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL)); 1239dbcead6SJeremy L Thompson 1249dbcead6SJeremy L Thompson PetscOptionsEnd(); 1259dbcead6SJeremy L Thompson 1269dbcead6SJeremy L Thompson // Create background mesh 1279dbcead6SJeremy L Thompson { 1289dbcead6SJeremy L Thompson PetscCall(DMCreate(comm, &dm_mesh)); 1299dbcead6SJeremy L Thompson PetscCall(DMSetType(dm_mesh, DMPLEX)); 1309dbcead6SJeremy L Thompson PetscCall(DMSetFromOptions(dm_mesh)); 1319dbcead6SJeremy L Thompson 1329dbcead6SJeremy L Thompson // -- Check for tensor product mesh 1339dbcead6SJeremy L Thompson { 1349dbcead6SJeremy L Thompson PetscBool is_simplex; 1359dbcead6SJeremy L Thompson 1369dbcead6SJeremy L Thompson PetscCall(DMPlexIsSimplex(dm_mesh, &is_simplex)); 1379dbcead6SJeremy L Thompson PetscCheck(!is_simplex, comm, PETSC_ERR_USER, "Only tensor-product background meshes supported"); 1389dbcead6SJeremy L Thompson } 1399dbcead6SJeremy L Thompson 1409dbcead6SJeremy L Thompson // -- Mesh FE space 1419dbcead6SJeremy L Thompson PetscCall(DMGetDimension(dm_mesh, &dim)); 1429dbcead6SJeremy L Thompson { 1439dbcead6SJeremy L Thompson PetscFE fe; 1449dbcead6SJeremy L Thompson 1459dbcead6SJeremy L Thompson PetscCall(DMGetDimension(dm_mesh, &dim)); 1469dbcead6SJeremy L Thompson PetscCall(PetscFECreateLagrange(comm, dim, num_comp, PETSC_FALSE, solution_order, solution_order + q_extra, &fe)); 1479dbcead6SJeremy L Thompson PetscCall(DMAddField(dm_mesh, NULL, (PetscObject)fe)); 1489dbcead6SJeremy L Thompson PetscCall(PetscFEDestroy(&fe)); 1499dbcead6SJeremy L Thompson } 1509dbcead6SJeremy L Thompson PetscCall(DMCreateDS(dm_mesh)); 1519dbcead6SJeremy L Thompson 1529dbcead6SJeremy L Thompson // -- Coordinate FE space 1539dbcead6SJeremy L Thompson { 1549dbcead6SJeremy L Thompson PetscFE fe_coord; 1559dbcead6SJeremy L Thompson 1569dbcead6SJeremy L Thompson PetscCall(PetscFECreateLagrange(comm, dim, dim, PETSC_FALSE, mesh_order, solution_order + q_extra, &fe_coord)); 15749a40c8aSKenneth E. Jansen PetscCall(DMSetCoordinateDisc(dm_mesh, fe_coord, PETSC_TRUE)); 1589dbcead6SJeremy L Thompson PetscCall(PetscFEDestroy(&fe_coord)); 1599dbcead6SJeremy L Thompson } 1609dbcead6SJeremy L Thompson 1619dbcead6SJeremy L Thompson // -- Set tensor permutation 1629dbcead6SJeremy L Thompson { 1639dbcead6SJeremy L Thompson DM dm_coord; 1649dbcead6SJeremy L Thompson 1659dbcead6SJeremy L Thompson PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord)); 1669dbcead6SJeremy L Thompson PetscCall(DMPlexSetClosurePermutationTensor(dm_mesh, PETSC_DETERMINE, NULL)); 1679dbcead6SJeremy L Thompson PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 1689dbcead6SJeremy L Thompson } 1699dbcead6SJeremy L Thompson 1709dbcead6SJeremy L Thompson // -- Final background mesh 1719dbcead6SJeremy L Thompson PetscCall(PetscObjectSetName((PetscObject)dm_mesh, "Background Mesh")); 1729dbcead6SJeremy L Thompson PetscCall(DMViewFromOptions(dm_mesh, NULL, "-dm_mesh_view")); 1739dbcead6SJeremy L Thompson } 1749dbcead6SJeremy L Thompson 1759dbcead6SJeremy L Thompson // Create particle swarm 1769dbcead6SJeremy L Thompson { 1779dbcead6SJeremy L Thompson PetscCall(DMCreate(comm, &dm_swarm)); 1789dbcead6SJeremy L Thompson PetscCall(DMSetType(dm_swarm, DMSWARM)); 1799dbcead6SJeremy L Thompson PetscCall(DMSetDimension(dm_swarm, dim)); 1809dbcead6SJeremy L Thompson PetscCall(DMSwarmSetType(dm_swarm, DMSWARM_PIC)); 1819dbcead6SJeremy L Thompson PetscCall(DMSwarmSetCellDM(dm_swarm, dm_mesh)); 1829dbcead6SJeremy L Thompson 1839dbcead6SJeremy L Thompson // -- Swarm field 1849dbcead6SJeremy L Thompson PetscCall(DMSwarmRegisterPetscDatatypeField(dm_swarm, DMSwarmPICField_u, num_comp, PETSC_SCALAR)); 1859dbcead6SJeremy L Thompson PetscCall(DMSwarmFinalizeFieldRegister(dm_swarm)); 1869dbcead6SJeremy L Thompson PetscCall(DMSwarmSetLocalSizes(dm_swarm, num_points, 0)); 1879dbcead6SJeremy L Thompson PetscCall(DMSetFromOptions(dm_swarm)); 1889dbcead6SJeremy L Thompson 1899dbcead6SJeremy L Thompson // -- Set swarm point locations 1909dbcead6SJeremy L Thompson PetscCall(DMSwarmInitalizePointLocations(dm_swarm, point_swarm_type, num_points, num_points_per_cell)); 1919dbcead6SJeremy L Thompson 1929dbcead6SJeremy L Thompson // -- Final particle swarm 1939dbcead6SJeremy L Thompson PetscCall(PetscObjectSetName((PetscObject)dm_swarm, "Particle Swarm")); 1949dbcead6SJeremy L Thompson PetscCall(DMViewFromOptions(dm_swarm, NULL, "-dm_swarm_view")); 1959dbcead6SJeremy L Thompson } 1969dbcead6SJeremy L Thompson 1979dbcead6SJeremy L Thompson // Set field values on background mesh 1989dbcead6SJeremy L Thompson PetscCall(DMCreateGlobalVector(dm_mesh, &U_mesh)); 1999dbcead6SJeremy L Thompson switch (target_type) { 2009dbcead6SJeremy L Thompson case TARGET_TANH: 2019dbcead6SJeremy L Thompson target_function_proj = EvalU_Tanh_proj; 2029dbcead6SJeremy L Thompson break; 2039dbcead6SJeremy L Thompson case TARGET_POLYNOMIAL: 2049dbcead6SJeremy L Thompson target_function_proj = EvalU_Poly_proj; 2059dbcead6SJeremy L Thompson break; 2069dbcead6SJeremy L Thompson case TARGET_SPHERE: 2079dbcead6SJeremy L Thompson target_function_proj = EvalU_Sphere_proj; 2089dbcead6SJeremy L Thompson break; 2099dbcead6SJeremy L Thompson } 2109dbcead6SJeremy L Thompson { 2119dbcead6SJeremy L Thompson TargetFuncProj mesh_solution[1] = {target_function_proj}; 2129dbcead6SJeremy L Thompson 2139dbcead6SJeremy L Thompson PetscCall(DMProjectFunction(dm_mesh, 0.0, mesh_solution, NULL, INSERT_VALUES, U_mesh)); 2149dbcead6SJeremy L Thompson } 2159dbcead6SJeremy L Thompson 2169dbcead6SJeremy L Thompson // Visualize background mesh 2179dbcead6SJeremy L Thompson PetscCall(PetscObjectSetName((PetscObject)U_mesh, "U on Background Mesh")); 2189dbcead6SJeremy L Thompson PetscCall(VecViewFromOptions(U_mesh, NULL, "-u_mesh_view")); 2199dbcead6SJeremy L Thompson 2209dbcead6SJeremy L Thompson // libCEED objects for swarm and background mesh 2219dbcead6SJeremy L Thompson PetscCall(DMSwarmCeedContextCreate(dm_swarm, ceed_resource, &swarm_ceed_context)); 2229dbcead6SJeremy L Thompson 2239dbcead6SJeremy L Thompson // Interpolate from mesh to points via PETSc 2249dbcead6SJeremy L Thompson { 2259dbcead6SJeremy L Thompson PetscCall(DMSwarmInterpolateFromCellToSwarm_Petsc(dm_swarm, DMSwarmPICField_u, U_mesh)); 2269dbcead6SJeremy L Thompson if (view_petsc_swarm) PetscCall(DMSwarmViewXDMF(dm_swarm, "swarm_petsc.xmf")); 2279dbcead6SJeremy L Thompson PetscCall(DMSwarmCheckSwarmValues(dm_swarm, DMSwarmPICField_u, tolerance, target_function_proj)); 2289dbcead6SJeremy L Thompson } 2299dbcead6SJeremy L Thompson 2309dbcead6SJeremy L Thompson // Interpolate from mesh to points via libCEED 2319dbcead6SJeremy L Thompson { 2329dbcead6SJeremy L Thompson PetscCall(DMSwarmInterpolateFromCellToSwarm_Ceed(dm_swarm, DMSwarmPICField_u, U_mesh)); 2339dbcead6SJeremy L Thompson if (view_ceed_swarm) PetscCall(DMSwarmViewXDMF(dm_swarm, "swarm_ceed.xmf")); 2349dbcead6SJeremy L Thompson PetscCall(DMSwarmCheckSwarmValues(dm_swarm, DMSwarmPICField_u, tolerance, target_function_proj)); 2359dbcead6SJeremy L Thompson } 2369dbcead6SJeremy L Thompson 2379dbcead6SJeremy L Thompson // Project from points to mesh via libCEED 2389dbcead6SJeremy L Thompson { 239*b6972d74SZach Atkins Vec U_projected; 2409dbcead6SJeremy L Thompson 2419dbcead6SJeremy L Thompson PetscCall(VecDuplicate(U_mesh, &U_projected)); 242*b6972d74SZach Atkins PetscCall(DMSwarmProjectFromSwarmToCells(dm_swarm, DMSwarmPICField_u, U_projected)); 2439dbcead6SJeremy L Thompson 244*b6972d74SZach Atkins PetscCall(PetscObjectSetName((PetscObject)U_projected, "U projected to Background Mesh")); 2459dbcead6SJeremy L Thompson PetscCall(VecViewFromOptions(U_projected, NULL, "-u_projected_view")); 2469dbcead6SJeremy L Thompson PetscCall(VecAXPY(U_projected, -1.0, U_mesh)); 2479dbcead6SJeremy L Thompson PetscCall(PetscObjectSetName((PetscObject)U_projected, "U Projection Error")); 2489dbcead6SJeremy L Thompson PetscCall(VecViewFromOptions(U_projected, NULL, "-u_error_view")); 249*b6972d74SZach Atkins 250*b6972d74SZach Atkins // -- Check error 2519dbcead6SJeremy L Thompson { 2529dbcead6SJeremy L Thompson PetscScalar error, norm_u_mesh; 2539dbcead6SJeremy L Thompson 2549dbcead6SJeremy L Thompson PetscCall(VecNorm(U_projected, NORM_2, &error)); 2559dbcead6SJeremy L Thompson PetscCall(VecNorm(U_mesh, NORM_2, &norm_u_mesh)); 2569dbcead6SJeremy L Thompson PetscCheck(error / norm_u_mesh < tolerance, comm, PETSC_ERR_USER, "Projection error too high: %e\n", error / norm_u_mesh); 2579dbcead6SJeremy L Thompson if (!test_mode) PetscCall(PetscPrintf(comm, " Projection error: %e\n", error / norm_u_mesh)); 2589dbcead6SJeremy L Thompson } 2599dbcead6SJeremy L Thompson 2609dbcead6SJeremy L Thompson PetscCall(VecDestroy(&U_projected)); 2619dbcead6SJeremy L Thompson } 2629dbcead6SJeremy L Thompson // Cleanup 2639dbcead6SJeremy L Thompson PetscCall(DMSwarmCeedContextDestroy(&swarm_ceed_context)); 2649dbcead6SJeremy L Thompson PetscCall(DMDestroy(&dm_swarm)); 2659dbcead6SJeremy L Thompson PetscCall(DMDestroy(&dm_mesh)); 2669dbcead6SJeremy L Thompson PetscCall(VecDestroy(&U_mesh)); 2679dbcead6SJeremy L Thompson return PetscFinalize(); 2689dbcead6SJeremy L Thompson } 2699dbcead6SJeremy L Thompson 2709dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 2719dbcead6SJeremy L Thompson // Solution functions 2729dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 2739dbcead6SJeremy L Thompson PetscScalar EvalU_Poly(PetscInt dim, const PetscScalar x[]) { 2749dbcead6SJeremy L Thompson PetscScalar result = 0.0; 2759dbcead6SJeremy L Thompson const PetscScalar p[5] = {3, 1, 4, 1, 5}; 2769dbcead6SJeremy L Thompson 2779dbcead6SJeremy L Thompson for (PetscInt d = 0; d < dim; d++) { 2789dbcead6SJeremy L Thompson PetscScalar result_1d = 1.0; 2799dbcead6SJeremy L Thompson 2809dbcead6SJeremy L Thompson for (PetscInt i = 4; i >= 0; i--) result_1d = result_1d * x[d] + p[i]; 2819dbcead6SJeremy L Thompson result += result_1d; 2829dbcead6SJeremy L Thompson } 2839dbcead6SJeremy L Thompson return result * 1E-3; 2849dbcead6SJeremy L Thompson } 2859dbcead6SJeremy L Thompson 2869dbcead6SJeremy L Thompson PetscScalar EvalU_Tanh(PetscInt dim, const PetscScalar x[]) { 2879dbcead6SJeremy L Thompson PetscScalar result = 1.0, center = 0.1; 2889dbcead6SJeremy L Thompson 2899dbcead6SJeremy L Thompson for (PetscInt d = 0; d < dim; d++) { 2909dbcead6SJeremy L Thompson result *= tanh(x[d] - center); 2919dbcead6SJeremy L Thompson center += 0.1; 2929dbcead6SJeremy L Thompson } 2939dbcead6SJeremy L Thompson return result; 2949dbcead6SJeremy L Thompson } 2959dbcead6SJeremy L Thompson 2969dbcead6SJeremy L Thompson PetscScalar EvalU_Sphere(PetscInt dim, const PetscScalar x[]) { 2979dbcead6SJeremy L Thompson PetscScalar distance = sqrt(x[0] * x[0] + x[1] * x[1] + x[2] * x[2]); 2989dbcead6SJeremy L Thompson 2999dbcead6SJeremy L Thompson return distance < 1.0 ? 1.0 : 0.1; 3009dbcead6SJeremy L Thompson } 3019dbcead6SJeremy L Thompson 3029dbcead6SJeremy L Thompson PetscErrorCode EvalU_Poly_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx) { 3039dbcead6SJeremy L Thompson PetscFunctionBeginUser; 3049dbcead6SJeremy L Thompson 3059dbcead6SJeremy L Thompson const PetscScalar f_x = EvalU_Poly(dim, x); 3069dbcead6SJeremy L Thompson 3079dbcead6SJeremy L Thompson for (PetscInt c = 0; c < num_comp; c++) u[c] = (c + 1.0) * f_x; 3089dbcead6SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 3099dbcead6SJeremy L Thompson } 3109dbcead6SJeremy L Thompson 3119dbcead6SJeremy L Thompson PetscErrorCode EvalU_Tanh_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx) { 3129dbcead6SJeremy L Thompson PetscFunctionBeginUser; 3139dbcead6SJeremy L Thompson 3149dbcead6SJeremy L Thompson const PetscScalar f_x = EvalU_Tanh(dim, x); 3159dbcead6SJeremy L Thompson 3169dbcead6SJeremy L Thompson for (PetscInt c = 0; c < num_comp; c++) u[c] = (c + 1.0) * f_x; 3179dbcead6SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 3189dbcead6SJeremy L Thompson } 3199dbcead6SJeremy L Thompson 3209dbcead6SJeremy L Thompson PetscErrorCode EvalU_Sphere_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx) { 3219dbcead6SJeremy L Thompson PetscFunctionBeginUser; 3229dbcead6SJeremy L Thompson 3239dbcead6SJeremy L Thompson const PetscScalar f_x = EvalU_Sphere(dim, x); 3249dbcead6SJeremy L Thompson 3259dbcead6SJeremy L Thompson for (PetscInt c = 0; c < num_comp; c++) u[c] = (c + 1.0) * f_x; 3269dbcead6SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 3279dbcead6SJeremy L Thompson } 3289dbcead6SJeremy L Thompson 3299dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 3309dbcead6SJeremy L Thompson // Projection via PETSc 3319dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 3329dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Petsc(DM dm_swarm, const char *field, Vec U_mesh) { 3339dbcead6SJeremy L Thompson PetscInt dim, num_comp, cell_start, cell_end; 3349dbcead6SJeremy L Thompson PetscScalar *u_points; 3359dbcead6SJeremy L Thompson const PetscScalar *coords_points; 3369dbcead6SJeremy L Thompson const PetscReal v0_ref[3] = {-1.0, -1.0, -1.0}; 3379dbcead6SJeremy L Thompson DM dm_mesh; 3389dbcead6SJeremy L Thompson PetscSection section_u_mesh_loc; 3399dbcead6SJeremy L Thompson PetscDS ds; 3409dbcead6SJeremy L Thompson PetscFE fe; 3419dbcead6SJeremy L Thompson PetscFEGeom fe_geometry; 3429dbcead6SJeremy L Thompson PetscQuadrature quadrature; 3439dbcead6SJeremy L Thompson Vec U_loc; 3449dbcead6SJeremy L Thompson 3459dbcead6SJeremy L Thompson PetscFunctionBeginUser; 3469dbcead6SJeremy L Thompson // Get mesh DM 3479dbcead6SJeremy L Thompson PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 3489dbcead6SJeremy L Thompson PetscCall(DMGetDimension(dm_mesh, &dim)); 3499dbcead6SJeremy L Thompson { 3509dbcead6SJeremy L Thompson PetscSection section_u_mesh_loc_closure_permutation; 3519dbcead6SJeremy L Thompson 3529dbcead6SJeremy L Thompson PetscCall(DMGetLocalSection(dm_mesh, §ion_u_mesh_loc_closure_permutation)); 3539dbcead6SJeremy L Thompson PetscCall(PetscSectionClone(section_u_mesh_loc_closure_permutation, §ion_u_mesh_loc)); 3549dbcead6SJeremy L Thompson PetscCall(PetscSectionResetClosurePermutation(section_u_mesh_loc)); 3559dbcead6SJeremy L Thompson } 3569dbcead6SJeremy L Thompson 3579dbcead6SJeremy L Thompson // Get local mesh values 3589dbcead6SJeremy L Thompson PetscCall(DMGetLocalVector(dm_mesh, &U_loc)); 3599dbcead6SJeremy L Thompson PetscCall(VecZeroEntries(U_loc)); 3609dbcead6SJeremy L Thompson PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_loc)); 3619dbcead6SJeremy L Thompson 3629dbcead6SJeremy L Thompson // Get local swarm data 3639dbcead6SJeremy L Thompson PetscCall(DMSwarmSortGetAccess(dm_swarm)); 3649dbcead6SJeremy L Thompson PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end)); 3659dbcead6SJeremy L Thompson PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points)); 3669dbcead6SJeremy L Thompson PetscCall(DMSwarmGetField(dm_swarm, field, &num_comp, NULL, (void **)&u_points)); 3679dbcead6SJeremy L Thompson 3689dbcead6SJeremy L Thompson // Interpolate values to each swarm point, one element in the background mesh at a time 3699dbcead6SJeremy L Thompson PetscCall(DMGetDS(dm_mesh, &ds)); 3709dbcead6SJeremy L Thompson PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe)); 3719dbcead6SJeremy L Thompson for (PetscInt cell = cell_start; cell < cell_end; cell++) { 3729dbcead6SJeremy L Thompson PetscTabulation tabulation; 3739dbcead6SJeremy L Thompson PetscScalar *u_cell = NULL, *coords_points_cell_true, *coords_points_cell_ref; 3749dbcead6SJeremy L Thompson PetscReal v[dim], J[dim * dim], invJ[dim * dim], detJ; 3759dbcead6SJeremy L Thompson PetscInt *points_cell; 3769dbcead6SJeremy L Thompson PetscInt num_points_in_cell; 3779dbcead6SJeremy L Thompson 3789dbcead6SJeremy L Thompson PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points_cell)); 3799dbcead6SJeremy L Thompson PetscCall(DMGetWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_true)); 3809dbcead6SJeremy L Thompson PetscCall(DMGetWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_ref)); 3819dbcead6SJeremy L Thompson // -- Reference coordinates for swarm points in background mesh element 3829dbcead6SJeremy L Thompson for (PetscInt p = 0; p < num_points_in_cell; p++) { 3839dbcead6SJeremy L Thompson for (PetscInt d = 0; d < dim; d++) coords_points_cell_true[p * dim + d] = coords_points[points_cell[p] * dim + d]; 3849dbcead6SJeremy L Thompson } 3859dbcead6SJeremy L Thompson PetscCall(DMPlexComputeCellGeometryFEM(dm_mesh, cell, NULL, v, J, invJ, &detJ)); 3869dbcead6SJeremy L Thompson for (PetscInt p = 0; p < num_points_in_cell; p++) { 3879dbcead6SJeremy L Thompson CoordinatesRealToRef(dim, dim, v0_ref, v, invJ, &coords_points_cell_true[p * dim], &coords_points_cell_ref[p * dim]); 3889dbcead6SJeremy L Thompson } 3899dbcead6SJeremy L Thompson // -- Interpolate values from current element in background mesh to swarm points 3909dbcead6SJeremy L Thompson PetscCall(PetscFECreateTabulation(fe, 1, num_points_in_cell, coords_points_cell_ref, 1, &tabulation)); 3919dbcead6SJeremy L Thompson PetscCall(DMPlexVecGetClosure(dm_mesh, section_u_mesh_loc, U_loc, cell, NULL, &u_cell)); 3929dbcead6SJeremy L Thompson PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 3939dbcead6SJeremy L Thompson PetscCall(PetscFECreateCellGeometry(fe, quadrature, &fe_geometry)); 3949dbcead6SJeremy L Thompson for (PetscInt p = 0; p < num_points_in_cell; p++) { 3959dbcead6SJeremy L Thompson PetscCall(PetscFEInterpolateAtPoints_Static(fe, tabulation, u_cell, &fe_geometry, p, &u_points[points_cell[p] * num_comp])); 3969dbcead6SJeremy L Thompson } 3979dbcead6SJeremy L Thompson 3989dbcead6SJeremy L Thompson // -- Cleanup 3999dbcead6SJeremy L Thompson PetscCall(PetscFEDestroyCellGeometry(fe, &fe_geometry)); 4009dbcead6SJeremy L Thompson PetscCall(DMPlexVecRestoreClosure(dm_mesh, section_u_mesh_loc, U_loc, cell, NULL, &u_cell)); 4019dbcead6SJeremy L Thompson PetscCall(DMRestoreWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_true)); 4029dbcead6SJeremy L Thompson PetscCall(DMRestoreWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_ref)); 4039dbcead6SJeremy L Thompson PetscCall(PetscTabulationDestroy(&tabulation)); 4049dbcead6SJeremy L Thompson PetscCall(PetscFree(points_cell)); 4059dbcead6SJeremy L Thompson } 4069dbcead6SJeremy L Thompson 4079dbcead6SJeremy L Thompson // Cleanup 4089dbcead6SJeremy L Thompson PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points)); 4099dbcead6SJeremy L Thompson PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&u_points)); 4109dbcead6SJeremy L Thompson PetscCall(DMSwarmSortRestoreAccess(dm_swarm)); 4119dbcead6SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm_mesh, &U_loc)); 4129dbcead6SJeremy L Thompson PetscCall(PetscSectionDestroy(§ion_u_mesh_loc)); 4139dbcead6SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 4149dbcead6SJeremy L Thompson } 4159dbcead6SJeremy L Thompson 4169dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 4179dbcead6SJeremy L Thompson // Projection via libCEED 4189dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 4199dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Ceed(DM dm_swarm, const char *field, Vec U_mesh) { 4209dbcead6SJeremy L Thompson PetscMemType U_mem_type; 4219dbcead6SJeremy L Thompson DM dm_mesh; 4229dbcead6SJeremy L Thompson Vec U_mesh_loc; 4239dbcead6SJeremy L Thompson DMSwarmCeedContext swarm_ceed_context; 4249dbcead6SJeremy L Thompson 4259dbcead6SJeremy L Thompson PetscFunctionBeginUser; 4269dbcead6SJeremy L Thompson // Get mesh DM 4279dbcead6SJeremy L Thompson PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 4289dbcead6SJeremy L Thompson PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context)); 4299dbcead6SJeremy L Thompson 4309dbcead6SJeremy L Thompson // Get mesh values 4319dbcead6SJeremy L Thompson PetscCall(DMGetLocalVector(dm_mesh, &U_mesh_loc)); 4329dbcead6SJeremy L Thompson PetscCall(VecZeroEntries(U_mesh_loc)); 4339dbcead6SJeremy L Thompson PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_mesh_loc)); 434*b6972d74SZach Atkins PetscCall(VecReadP2C(U_mesh_loc, &U_mem_type, swarm_ceed_context->u_mesh)); 4359dbcead6SJeremy L Thompson 4369dbcead6SJeremy L Thompson // Get swarm access 4379dbcead6SJeremy L Thompson PetscCall(DMSwarmSortGetAccess(dm_swarm)); 438*b6972d74SZach Atkins PetscCall(DMSwarmPICFieldP2C(dm_swarm, field, swarm_ceed_context->u_points)); 4399dbcead6SJeremy L Thompson 440*b6972d74SZach Atkins // Interpolate field from mesh to swarm points 441*b6972d74SZach Atkins CeedOperatorApply(swarm_ceed_context->op_mesh_to_points, swarm_ceed_context->u_mesh, swarm_ceed_context->u_points, CEED_REQUEST_IMMEDIATE); 4429dbcead6SJeremy L Thompson 4439dbcead6SJeremy L Thompson // Cleanup 444*b6972d74SZach Atkins PetscCall(DMSwarmPICFieldC2P(dm_swarm, field, swarm_ceed_context->u_points)); 4459dbcead6SJeremy L Thompson PetscCall(DMSwarmSortRestoreAccess(dm_swarm)); 446*b6972d74SZach Atkins PetscCall(VecReadC2P(swarm_ceed_context->u_mesh, U_mem_type, U_mesh_loc)); 4479dbcead6SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm_mesh, &U_mesh_loc)); 4489dbcead6SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 4499dbcead6SJeremy L Thompson } 4509dbcead6SJeremy L Thompson 4519dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 4529dbcead6SJeremy L Thompson // Error checking utility 4539dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------ 4549dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCheckSwarmValues(DM dm_swarm, const char *field, PetscScalar tolerance, TargetFuncProj TrueSolution) { 4559dbcead6SJeremy L Thompson PetscBool within_tolerance = PETSC_TRUE; 4569dbcead6SJeremy L Thompson PetscInt dim, num_comp, cell_start, cell_end; 4579dbcead6SJeremy L Thompson const PetscScalar *u_points, *coords_points; 4589dbcead6SJeremy L Thompson DM dm_mesh; 4599dbcead6SJeremy L Thompson 4609dbcead6SJeremy L Thompson PetscFunctionBeginUser; 4619dbcead6SJeremy L Thompson PetscCall(DMSwarmSortGetAccess(dm_swarm)); 4629dbcead6SJeremy L Thompson PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 4639dbcead6SJeremy L Thompson PetscCall(DMGetDimension(dm_mesh, &dim)); 4649dbcead6SJeremy L Thompson PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end)); 4659dbcead6SJeremy L Thompson PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points)); 4669dbcead6SJeremy L Thompson PetscCall(DMSwarmGetField(dm_swarm, field, &num_comp, NULL, (void **)&u_points)); 4679dbcead6SJeremy L Thompson 4689dbcead6SJeremy L Thompson // Interpolate values to each swarm point, one element in the background mesh at a time 4699dbcead6SJeremy L Thompson for (PetscInt cell = cell_start; cell < cell_end; cell++) { 4709dbcead6SJeremy L Thompson PetscInt *points; 4719dbcead6SJeremy L Thompson PetscInt num_points_in_cell; 4729dbcead6SJeremy L Thompson 4739dbcead6SJeremy L Thompson PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points)); 4749dbcead6SJeremy L Thompson // -- Reference coordinates for swarm points in background mesh element 4759dbcead6SJeremy L Thompson for (PetscInt p = 0; p < num_points_in_cell; p++) { 4769dbcead6SJeremy L Thompson PetscScalar x[dim], u_true[num_comp]; 4779dbcead6SJeremy L Thompson 4789dbcead6SJeremy L Thompson for (PetscInt d = 0; d < dim; d++) x[d] = coords_points[points[p] * dim + d]; 4799dbcead6SJeremy L Thompson PetscCall(TrueSolution(dim, 0.0, x, num_comp, u_true, NULL)); 4809dbcead6SJeremy L Thompson for (PetscInt i = 0; i < num_comp; i++) { 4819dbcead6SJeremy L Thompson if (PetscAbs(u_points[points[p] * num_comp + i] - u_true[i]) > tolerance) { 4829dbcead6SJeremy L Thompson within_tolerance = PETSC_FALSE; 4839dbcead6SJeremy L Thompson PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm_swarm), 4849dbcead6SJeremy L Thompson "Incorrect interpolated value, cell %" PetscInt_FMT " point %" PetscInt_FMT " component %" PetscInt_FMT 4859dbcead6SJeremy L Thompson ", found %f expected %f\n", 4869dbcead6SJeremy L Thompson cell, p, i, u_points[points[p] * num_comp + i], u_true[i])); 4879dbcead6SJeremy L Thompson } 4889dbcead6SJeremy L Thompson } 4899dbcead6SJeremy L Thompson } 4909dbcead6SJeremy L Thompson 4919dbcead6SJeremy L Thompson // -- Cleanup 4929dbcead6SJeremy L Thompson PetscCall(PetscFree(points)); 4939dbcead6SJeremy L Thompson } 4949dbcead6SJeremy L Thompson 4959dbcead6SJeremy L Thompson // Cleanup 4969dbcead6SJeremy L Thompson PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points)); 4979dbcead6SJeremy L Thompson PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&u_points)); 4989dbcead6SJeremy L Thompson PetscCall(DMSwarmSortRestoreAccess(dm_swarm)); 4999dbcead6SJeremy L Thompson PetscCheck(within_tolerance, PetscObjectComm((PetscObject)dm_swarm), PETSC_ERR_USER, "Interpolation to swarm points not within tolerance"); 5009dbcead6SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 5019dbcead6SJeremy L Thompson } 502