xref: /libCEED/examples/petsc/include/swarmutils.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2b6972d74SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3b6972d74SZach Atkins //
4b6972d74SZach Atkins // SPDX-License-Identifier: BSD-2-Clause
5b6972d74SZach Atkins //
6b6972d74SZach Atkins // This file is part of CEED:  http://github.com/ceed
7b6972d74SZach Atkins 
8b6972d74SZach Atkins /// @file
9b6972d74SZach Atkins /// Utility functions for particle-based methods with DMSwarm
10b6972d74SZach Atkins #pragma once
11b6972d74SZach Atkins 
12b6972d74SZach Atkins #include <ceed.h>
13b6972d74SZach Atkins #include <math.h>
14b6972d74SZach Atkins #include <petscdmplex.h>
15b6972d74SZach Atkins #include <petscdmswarm.h>
16b6972d74SZach Atkins #include <petsc/private/petscfeimpl.h> /* For interpolation */
17b6972d74SZach Atkins 
18b6972d74SZach Atkins #include "petscutils.h"
19b6972d74SZach Atkins 
20b6972d74SZach Atkins // libCEED context data
21b6972d74SZach Atkins typedef struct DMSwarmCeedContext_ *DMSwarmCeedContext;
22b6972d74SZach Atkins struct DMSwarmCeedContext_ {
23b6972d74SZach Atkins   Ceed         ceed;
24b6972d74SZach Atkins   CeedVector   u_mesh, v_mesh, u_points;
25b6972d74SZach Atkins   CeedOperator op_mass, op_mesh_to_points, op_points_to_mesh;
26b6972d74SZach Atkins };
27b6972d74SZach Atkins 
28b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextCreate(DM dm_swarm, const char *ceed_resource, DMSwarmCeedContext *ctx);
29b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextDestroy(DMSwarmCeedContext *ctx);
30b6972d74SZach Atkins 
31b6972d74SZach Atkins // Swarm point distribution
32b6972d74SZach Atkins typedef enum { SWARM_GAUSS = 0, SWARM_UNIFORM = 1, SWARM_CELL_RANDOM = 2, SWARM_SINUSOIDAL = 3 } PointSwarmType;
33b6972d74SZach Atkins static const char *const point_swarm_types[] = {"gauss", "uniform", "cell_random", "sinusoidal", "PointSwarmType", "SWARM", 0};
34b6972d74SZach Atkins 
35b6972d74SZach Atkins // Memory utilities
36b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldP2C(DM dm_swarm, const char *field, CeedVector x_ceed);
37b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldC2P(DM dm_swarm, const char *field, CeedVector x_ceed);
38b6972d74SZach Atkins 
39b6972d74SZach Atkins // Swarm helper function
40b6972d74SZach Atkins PetscErrorCode DMSwarmInitalizePointLocations(DM dm_swarm, PointSwarmType point_swarm_type, PetscInt num_points, PetscInt num_points_per_cell);
41b6972d74SZach Atkins PetscErrorCode DMSwarmCreateReferenceCoordinates(DM dm_swarm, IS *is_points, Vec *ref_coords);
42b6972d74SZach Atkins 
43b6972d74SZach Atkins // Swarm to mesh projection
4478f7fce3SZach Atkins PetscErrorCode DMSwarmCreateProjectionRHS(DM dm_swarm, const char *field, Vec U_points, Vec B_mesh);
45b6972d74SZach Atkins PetscErrorCode MatMult_SwarmMass(Mat A, Vec U_mesh, Vec V_mesh);
4678f7fce3SZach Atkins PetscErrorCode DMSwarmProjectFromSwarmToCells(DM dm_swarm, const char *field, Vec U_points, Vec U_mesh);
4778f7fce3SZach Atkins 
4878f7fce3SZach Atkins PetscErrorCode SetupProblemSwarm(DM dm_swarm, Ceed ceed, BPData bp_data, CeedData data, PetscBool setup_rhs, Vec rhs, Vec target);
49