xref: /honee/src/bc_definition.c (revision 09240e0a7b1ed022a0b5d6f639eacef49a47bbbd)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3487a3b6eSJames Wright 
4487a3b6eSJames Wright #include <bc_definition.h>
5*09240e0aSJames Wright #include <petsc/private/petscimpl.h>
6487a3b6eSJames Wright 
7487a3b6eSJames Wright /**
8487a3b6eSJames Wright    @brief Create `BCDefinition`
9487a3b6eSJames Wright 
10487a3b6eSJames Wright    @param[in]  name             Name of the boundary condition
11487a3b6eSJames Wright    @param[in]  num_label_values Number of `DMLabel` values
12487a3b6eSJames Wright    @param[in]  label_values     Array of label values that define the boundaries controlled by the `BCDefinition`, size `num_label_values`
13487a3b6eSJames Wright    @param[out] bc_def           The new `BCDefinition`
14487a3b6eSJames Wright **/
15487a3b6eSJames Wright PetscErrorCode BCDefinitionCreate(const char *name, PetscInt num_label_values, PetscInt label_values[], BCDefinition *bc_def) {
16487a3b6eSJames Wright   PetscFunctionBeginUser;
17487a3b6eSJames Wright   PetscCall(PetscNew(bc_def));
18487a3b6eSJames Wright 
19487a3b6eSJames Wright   PetscCall(PetscStrallocpy(name, &(*bc_def)->name));
20487a3b6eSJames Wright   (*bc_def)->num_label_values = num_label_values;
21487a3b6eSJames Wright   PetscCall(PetscMalloc1(num_label_values, &(*bc_def)->label_values));
22487a3b6eSJames Wright   for (PetscInt i = 0; i < num_label_values; i++) (*bc_def)->label_values[i] = label_values[i];
23487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
24487a3b6eSJames Wright }
25487a3b6eSJames Wright 
26487a3b6eSJames Wright /**
27487a3b6eSJames Wright    @brief Get base information for `BCDefinition`
28487a3b6eSJames Wright 
29487a3b6eSJames Wright    @param[in]  bc_def           `BCDefinition` to get information from
30487a3b6eSJames Wright    @param[out] name             Name of the `BCDefinition`
31487a3b6eSJames Wright    @param[out] num_label_values Number of `DMLabel` values
32487a3b6eSJames Wright    @param[out] label_values     Array of label values that define the boundaries controlled by the `BCDefinition`, size `num_label_values`
33487a3b6eSJames Wright **/
34487a3b6eSJames Wright PetscErrorCode BCDefinitionGetInfo(BCDefinition bc_def, const char *name[], PetscInt *num_label_values, const PetscInt *label_values[]) {
35487a3b6eSJames Wright   PetscFunctionBeginUser;
36487a3b6eSJames Wright   if (name) *name = bc_def->name;
37487a3b6eSJames Wright   if (label_values) {
38487a3b6eSJames Wright     *num_label_values = bc_def->num_label_values;
39487a3b6eSJames Wright     *label_values     = bc_def->label_values;
40487a3b6eSJames Wright   }
41487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
42487a3b6eSJames Wright }
43487a3b6eSJames Wright 
44487a3b6eSJames Wright /**
45487a3b6eSJames Wright    @brief Destory a `BCDefinition` object
46487a3b6eSJames Wright 
47487a3b6eSJames Wright    @param[in,out] bc_def `BCDefinition` to be destroyed
48487a3b6eSJames Wright **/
49487a3b6eSJames Wright PetscErrorCode BCDefinitionDestroy(BCDefinition *bc_def) {
50487a3b6eSJames Wright   PetscFunctionBeginUser;
51487a3b6eSJames Wright   if ((*bc_def)->name) PetscCall(PetscFree((*bc_def)->name));
52487a3b6eSJames Wright   if ((*bc_def)->label_values) PetscCall(PetscFree((*bc_def)->label_values));
53487a3b6eSJames Wright   if ((*bc_def)->essential_comps) PetscCall(PetscFree((*bc_def)->essential_comps));
54*09240e0aSJames Wright   if ((*bc_def)->dm) PetscCall(DMDestroy(&(*bc_def)->dm));
55487a3b6eSJames Wright   PetscCall(PetscFree(*bc_def));
56487a3b6eSJames Wright   *bc_def = NULL;
57487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
58487a3b6eSJames Wright }
59487a3b6eSJames Wright 
60487a3b6eSJames Wright /**
61487a3b6eSJames Wright    @brief Set `DM_BC_ESSENTIAL` boundary condition values
62487a3b6eSJames Wright 
63487a3b6eSJames Wright    @param[in,out] bc_def              `BCDefinition` to set values to
64487a3b6eSJames Wright    @param[in]     num_essential_comps Number of components to set
65487a3b6eSJames Wright    @param[in]     essential_comps     Array of components to set, size `num_essential_comps`
66487a3b6eSJames Wright **/
67487a3b6eSJames Wright PetscErrorCode BCDefinitionSetEssential(BCDefinition bc_def, PetscInt num_essential_comps, PetscInt essential_comps[]) {
68487a3b6eSJames Wright   PetscFunctionBeginUser;
69487a3b6eSJames Wright   bc_def->num_essential_comps = num_essential_comps;
70487a3b6eSJames Wright   PetscCall(PetscMalloc1(num_essential_comps, &bc_def->essential_comps));
71487a3b6eSJames Wright   PetscCall(PetscArraycpy(bc_def->essential_comps, essential_comps, num_essential_comps));
72487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
73487a3b6eSJames Wright }
74487a3b6eSJames Wright 
75487a3b6eSJames Wright /**
76487a3b6eSJames Wright    @brief Get `DM_BC_ESSENTIAL` boundary condition values
77487a3b6eSJames Wright 
78487a3b6eSJames Wright    @param[in]  bc_def              `BCDefinition` to set values to
79487a3b6eSJames Wright    @param[out] num_essential_comps Number of components to set
80487a3b6eSJames Wright    @param[out] essential_comps     Array of components to set, size `num_essential_comps`
81487a3b6eSJames Wright **/
82487a3b6eSJames Wright PetscErrorCode BCDefinitionGetEssential(BCDefinition bc_def, PetscInt *num_essential_comps, const PetscInt *essential_comps[]) {
83487a3b6eSJames Wright   PetscFunctionBeginUser;
84487a3b6eSJames Wright   *num_essential_comps = bc_def->num_essential_comps;
85487a3b6eSJames Wright   *essential_comps     = bc_def->essential_comps;
86487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
87487a3b6eSJames Wright }
88487a3b6eSJames Wright 
89487a3b6eSJames Wright #define LABEL_ARRAY_SIZE 256
90487a3b6eSJames Wright 
91487a3b6eSJames Wright // @brief See `PetscOptionsBCDefinition`
92ddf6e248SJames Wright PetscErrorCode PetscOptionsBCDefinition_Private(PetscOptionItems PetscOptionsObject, const char opt[], const char text[], const char man[],
93487a3b6eSJames Wright                                                 const char name[], BCDefinition *bc_def, PetscBool *set) {
94487a3b6eSJames Wright   PetscInt num_label_values = LABEL_ARRAY_SIZE, label_values[LABEL_ARRAY_SIZE] = {0};
95487a3b6eSJames Wright 
96487a3b6eSJames Wright   PetscFunctionBeginUser;
97487a3b6eSJames Wright   PetscCall(PetscOptionsIntArray(opt, text, man, label_values, &num_label_values, set));
98487a3b6eSJames Wright   if (num_label_values > 0) {
99487a3b6eSJames Wright     PetscCall(BCDefinitionCreate(name, num_label_values, label_values, bc_def));
100487a3b6eSJames Wright   } else {
101487a3b6eSJames Wright     *bc_def = NULL;
102487a3b6eSJames Wright   }
103487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
104487a3b6eSJames Wright }
105*09240e0aSJames Wright 
106*09240e0aSJames Wright /**
107*09240e0aSJames Wright   @brief Set `DM` for BCDefinition
108*09240e0aSJames Wright 
109*09240e0aSJames Wright   @param[in,out] bc_def `BCDefinition` to add `dm` to
110*09240e0aSJames Wright   @param[in]     dm     `DM` to assign to BCDefinition, or `NULL` to remove `DM`
111*09240e0aSJames Wright **/
112*09240e0aSJames Wright PetscErrorCode BCDefinitionSetDM(BCDefinition bc_def, DM dm) {
113*09240e0aSJames Wright   PetscFunctionBeginUser;
114*09240e0aSJames Wright   if (bc_def->dm) PetscCall(DMDestroy(&bc_def->dm));
115*09240e0aSJames Wright   if (dm) {
116*09240e0aSJames Wright     PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
117*09240e0aSJames Wright     PetscCall(PetscObjectReference((PetscObject)dm));
118*09240e0aSJames Wright     bc_def->dm = dm;
119*09240e0aSJames Wright   }
120*09240e0aSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
121*09240e0aSJames Wright }
122*09240e0aSJames Wright 
123*09240e0aSJames Wright /**
124*09240e0aSJames Wright   @brief Get `DM` assigned to BCDefinition
125*09240e0aSJames Wright 
126*09240e0aSJames Wright   @param[in]  bc_def `BCDefinition` to get `dm` from
127*09240e0aSJames Wright   @param[out] dm     `DM` assigned to BCDefinition
128*09240e0aSJames Wright **/
129*09240e0aSJames Wright PetscErrorCode BCDefinitionGetDM(BCDefinition bc_def, DM *dm) {
130*09240e0aSJames Wright   PetscFunctionBeginUser;
131*09240e0aSJames Wright   PetscAssertPointer(dm, 2);
132*09240e0aSJames Wright   *dm = bc_def->dm;
133*09240e0aSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
134*09240e0aSJames Wright }
135