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