xref: /honee/src/bc_definition.c (revision 77aa5ad7c9c7ea06edea746b0c60eea6fdad3cd2)
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*77aa5ad7SJames Wright #include <dm-utils.h>
6*77aa5ad7SJames Wright #include <petsc-ceed.h>
709240e0aSJames Wright #include <petsc/private/petscimpl.h>
8487a3b6eSJames Wright 
9487a3b6eSJames Wright /**
10487a3b6eSJames Wright    @brief Create `BCDefinition`
11487a3b6eSJames Wright 
12487a3b6eSJames Wright    @param[in]  name             Name of the boundary condition
13487a3b6eSJames Wright    @param[in]  num_label_values Number of `DMLabel` values
14487a3b6eSJames Wright    @param[in]  label_values     Array of label values that define the boundaries controlled by the `BCDefinition`, size `num_label_values`
15487a3b6eSJames Wright    @param[out] bc_def           The new `BCDefinition`
16487a3b6eSJames Wright **/
17487a3b6eSJames Wright PetscErrorCode BCDefinitionCreate(const char *name, PetscInt num_label_values, PetscInt label_values[], BCDefinition *bc_def) {
18487a3b6eSJames Wright   PetscFunctionBeginUser;
19487a3b6eSJames Wright   PetscCall(PetscNew(bc_def));
20487a3b6eSJames Wright 
21487a3b6eSJames Wright   PetscCall(PetscStrallocpy(name, &(*bc_def)->name));
22487a3b6eSJames Wright   (*bc_def)->num_label_values = num_label_values;
23487a3b6eSJames Wright   PetscCall(PetscMalloc1(num_label_values, &(*bc_def)->label_values));
24487a3b6eSJames Wright   for (PetscInt i = 0; i < num_label_values; i++) (*bc_def)->label_values[i] = label_values[i];
25487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
26487a3b6eSJames Wright }
27487a3b6eSJames Wright 
28487a3b6eSJames Wright /**
29487a3b6eSJames Wright    @brief Get base information for `BCDefinition`
30487a3b6eSJames Wright 
31487a3b6eSJames Wright    @param[in]  bc_def           `BCDefinition` to get information from
32487a3b6eSJames Wright    @param[out] name             Name of the `BCDefinition`
33487a3b6eSJames Wright    @param[out] num_label_values Number of `DMLabel` values
34487a3b6eSJames Wright    @param[out] label_values     Array of label values that define the boundaries controlled by the `BCDefinition`, size `num_label_values`
35487a3b6eSJames Wright **/
36487a3b6eSJames Wright PetscErrorCode BCDefinitionGetInfo(BCDefinition bc_def, const char *name[], PetscInt *num_label_values, const PetscInt *label_values[]) {
37487a3b6eSJames Wright   PetscFunctionBeginUser;
38487a3b6eSJames Wright   if (name) *name = bc_def->name;
39487a3b6eSJames Wright   if (label_values) {
40487a3b6eSJames Wright     *num_label_values = bc_def->num_label_values;
41487a3b6eSJames Wright     *label_values     = bc_def->label_values;
42487a3b6eSJames Wright   }
43487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
44487a3b6eSJames Wright }
45487a3b6eSJames Wright 
46487a3b6eSJames Wright /**
47487a3b6eSJames Wright    @brief Destory a `BCDefinition` object
48487a3b6eSJames Wright 
49487a3b6eSJames Wright    @param[in,out] bc_def `BCDefinition` to be destroyed
50487a3b6eSJames Wright **/
51487a3b6eSJames Wright PetscErrorCode BCDefinitionDestroy(BCDefinition *bc_def) {
522e678684SJames Wright   BCDefinition bc_def_ = *bc_def;
53487a3b6eSJames Wright   PetscFunctionBeginUser;
542e678684SJames Wright   if (bc_def_->name) PetscCall(PetscFree(bc_def_->name));
552e678684SJames Wright   if (bc_def_->label_values) PetscCall(PetscFree(bc_def_->label_values));
562e678684SJames Wright   if (bc_def_->essential_comps) PetscCall(PetscFree(bc_def_->essential_comps));
572e678684SJames Wright   if (bc_def_->dm) PetscCall(DMDestroy(&bc_def_->dm));
582e678684SJames Wright   if (bc_def_->DestroyCtx) PetscCall((*(bc_def_->DestroyCtx))(&bc_def_->ctx));
592e678684SJames Wright   PetscCall(PetscFree(bc_def_));
60487a3b6eSJames Wright   *bc_def = NULL;
61487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
62487a3b6eSJames Wright }
63487a3b6eSJames Wright 
64487a3b6eSJames Wright /**
65487a3b6eSJames Wright    @brief Set `DM_BC_ESSENTIAL` boundary condition values
66487a3b6eSJames Wright 
67487a3b6eSJames Wright    @param[in,out] bc_def              `BCDefinition` to set values to
68487a3b6eSJames Wright    @param[in]     num_essential_comps Number of components to set
69487a3b6eSJames Wright    @param[in]     essential_comps     Array of components to set, size `num_essential_comps`
70487a3b6eSJames Wright **/
71487a3b6eSJames Wright PetscErrorCode BCDefinitionSetEssential(BCDefinition bc_def, PetscInt num_essential_comps, PetscInt essential_comps[]) {
72487a3b6eSJames Wright   PetscFunctionBeginUser;
73487a3b6eSJames Wright   bc_def->num_essential_comps = num_essential_comps;
74487a3b6eSJames Wright   PetscCall(PetscMalloc1(num_essential_comps, &bc_def->essential_comps));
75487a3b6eSJames Wright   PetscCall(PetscArraycpy(bc_def->essential_comps, essential_comps, num_essential_comps));
76487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
77487a3b6eSJames Wright }
78487a3b6eSJames Wright 
79487a3b6eSJames Wright /**
80487a3b6eSJames Wright    @brief Get `DM_BC_ESSENTIAL` boundary condition values
81487a3b6eSJames Wright 
82487a3b6eSJames Wright    @param[in]  bc_def              `BCDefinition` to set values to
83487a3b6eSJames Wright    @param[out] num_essential_comps Number of components to set
84487a3b6eSJames Wright    @param[out] essential_comps     Array of components to set, size `num_essential_comps`
85487a3b6eSJames Wright **/
86487a3b6eSJames Wright PetscErrorCode BCDefinitionGetEssential(BCDefinition bc_def, PetscInt *num_essential_comps, const PetscInt *essential_comps[]) {
87487a3b6eSJames Wright   PetscFunctionBeginUser;
88487a3b6eSJames Wright   *num_essential_comps = bc_def->num_essential_comps;
89487a3b6eSJames Wright   *essential_comps     = bc_def->essential_comps;
90487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
91487a3b6eSJames Wright }
92487a3b6eSJames Wright 
93487a3b6eSJames Wright #define LABEL_ARRAY_SIZE 256
94487a3b6eSJames Wright 
95487a3b6eSJames Wright // @brief See `PetscOptionsBCDefinition`
96ddf6e248SJames Wright PetscErrorCode PetscOptionsBCDefinition_Private(PetscOptionItems PetscOptionsObject, const char opt[], const char text[], const char man[],
97487a3b6eSJames Wright                                                 const char name[], BCDefinition *bc_def, PetscBool *set) {
98487a3b6eSJames Wright   PetscInt num_label_values = LABEL_ARRAY_SIZE, label_values[LABEL_ARRAY_SIZE] = {0};
99487a3b6eSJames Wright 
100487a3b6eSJames Wright   PetscFunctionBeginUser;
101487a3b6eSJames Wright   PetscCall(PetscOptionsIntArray(opt, text, man, label_values, &num_label_values, set));
102487a3b6eSJames Wright   if (num_label_values > 0) {
103487a3b6eSJames Wright     PetscCall(BCDefinitionCreate(name, num_label_values, label_values, bc_def));
104487a3b6eSJames Wright   } else {
105487a3b6eSJames Wright     *bc_def = NULL;
106487a3b6eSJames Wright   }
107487a3b6eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
108487a3b6eSJames Wright }
10909240e0aSJames Wright 
11009240e0aSJames Wright /**
11109240e0aSJames Wright   @brief Set `DM` for BCDefinition
11209240e0aSJames Wright 
11309240e0aSJames Wright   @param[in,out] bc_def `BCDefinition` to add `dm` to
11409240e0aSJames Wright   @param[in]     dm     `DM` to assign to BCDefinition, or `NULL` to remove `DM`
11509240e0aSJames Wright **/
11609240e0aSJames Wright PetscErrorCode BCDefinitionSetDM(BCDefinition bc_def, DM dm) {
11709240e0aSJames Wright   PetscFunctionBeginUser;
11809240e0aSJames Wright   if (bc_def->dm) PetscCall(DMDestroy(&bc_def->dm));
11909240e0aSJames Wright   if (dm) {
12009240e0aSJames Wright     PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
12109240e0aSJames Wright     PetscCall(PetscObjectReference((PetscObject)dm));
12209240e0aSJames Wright     bc_def->dm = dm;
12309240e0aSJames Wright   }
12409240e0aSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
12509240e0aSJames Wright }
12609240e0aSJames Wright 
12709240e0aSJames Wright /**
12809240e0aSJames Wright   @brief Get `DM` assigned to BCDefinition
12909240e0aSJames Wright 
13009240e0aSJames Wright   @param[in]  bc_def `BCDefinition` to get `dm` from
13109240e0aSJames Wright   @param[out] dm     `DM` assigned to BCDefinition
13209240e0aSJames Wright **/
13309240e0aSJames Wright PetscErrorCode BCDefinitionGetDM(BCDefinition bc_def, DM *dm) {
13409240e0aSJames Wright   PetscFunctionBeginUser;
13509240e0aSJames Wright   PetscAssertPointer(dm, 2);
13609240e0aSJames Wright   *dm = bc_def->dm;
13709240e0aSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
13809240e0aSJames Wright }
1392e678684SJames Wright 
1402e678684SJames Wright /**
1412e678684SJames Wright   @brief Set custom context struct for use in BCDefinition
1422e678684SJames Wright 
1432e678684SJames Wright   @param[in,out] bc_def      `BCDefinition` to add `ctx` to
1442e678684SJames Wright   @param[in]     destroy_ctx Optional function pointer that destroys the user context on `BCDefinitionDestroy()`
1452e678684SJames Wright   @param[in]     ctx         Pointer to context struct
1462e678684SJames Wright **/
1472e678684SJames Wright PetscErrorCode BCDefinitionSetContext(BCDefinition bc_def, PetscCtxDestroyFn *destroy_ctx, void *ctx) {
1482e678684SJames Wright   PetscFunctionBeginUser;
1492e678684SJames Wright   if (bc_def->DestroyCtx) PetscCall((*(bc_def->DestroyCtx))(&bc_def->ctx));
1502e678684SJames Wright   bc_def->ctx        = ctx;
1512e678684SJames Wright   bc_def->DestroyCtx = destroy_ctx;
1522e678684SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1532e678684SJames Wright }
1542e678684SJames Wright 
1552e678684SJames Wright /**
1562e678684SJames Wright   @brief Set custom context struct for use in BCDefinition
1572e678684SJames Wright 
1582e678684SJames Wright   @param[in]  bc_def `BCDefinition` to get `ctx` from
1592e678684SJames Wright   @param[out] ctx    Pointer to context struct
1602e678684SJames Wright **/
1612e678684SJames Wright PetscErrorCode BCDefinitionGetContext(BCDefinition bc_def, void *ctx) {
1622e678684SJames Wright   PetscFunctionBeginUser;
1632e678684SJames Wright   PetscAssertPointer(ctx, 2);
1642e678684SJames Wright   *(void **)ctx = bc_def->ctx;
1652e678684SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1662e678684SJames Wright }
167*77aa5ad7SJames Wright 
168*77aa5ad7SJames Wright /**
169*77aa5ad7SJames Wright   @brief Add function pointers to create `CeedQFunction` and `CeedOperator` for IFunction of boundary condition
170*77aa5ad7SJames Wright 
171*77aa5ad7SJames Wright   @param[in,out] bc_def    `BCDefinition` to add function pointers to
172*77aa5ad7SJames Wright   @param[in]     create_qf Function to create `CeedQFunction`
173*77aa5ad7SJames Wright   @param[in]     add_op    Function to create and add `CeedOperator` to composite `CeedOperator`
174*77aa5ad7SJames Wright **/
175*77aa5ad7SJames Wright PetscErrorCode BCDefinitionSetIFunction(BCDefinition bc_def, BCDefinitionCreateQFunction create_qf, BCDefinitionAddIFunctionOperator add_op) {
176*77aa5ad7SJames Wright   PetscFunctionBeginUser;
177*77aa5ad7SJames Wright   bc_def->CreateIFunctionQF    = create_qf;
178*77aa5ad7SJames Wright   bc_def->AddIFunctionOperator = add_op;
179*77aa5ad7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
180*77aa5ad7SJames Wright }
181*77aa5ad7SJames Wright 
182*77aa5ad7SJames Wright /**
183*77aa5ad7SJames Wright   @brief Add function pointers to create `CeedQFunction` and `CeedOperator` for IJacobian of boundary condition
184*77aa5ad7SJames Wright 
185*77aa5ad7SJames Wright   @param[in,out] bc_def    `BCDefinition` to add function pointers to
186*77aa5ad7SJames Wright   @param[in]     create_qf Function to create `CeedQFunction`
187*77aa5ad7SJames Wright   @param[in]     add_op    Function to create and add `CeedOperator` to composite `CeedOperator`
188*77aa5ad7SJames Wright **/
189*77aa5ad7SJames Wright PetscErrorCode BCDefinitionSetIJacobian(BCDefinition bc_def, BCDefinitionCreateQFunction create_qf, BCDefinitionAddIJacobianOperator add_op) {
190*77aa5ad7SJames Wright   PetscFunctionBeginUser;
191*77aa5ad7SJames Wright   bc_def->CreateIJacobianQF    = create_qf;
192*77aa5ad7SJames Wright   bc_def->AddIJacobianOperator = add_op;
193*77aa5ad7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
194*77aa5ad7SJames Wright }
195*77aa5ad7SJames Wright 
196*77aa5ad7SJames Wright /**
197*77aa5ad7SJames Wright   @brief Add operators (IFunction, IJacobian) to composite operator
198*77aa5ad7SJames Wright 
199*77aa5ad7SJames Wright   This loops over orientations for each face label specified by `bc_def` and adds the IFunction and IJacobian operator to respective composite operator.
200*77aa5ad7SJames Wright 
201*77aa5ad7SJames Wright   @param[in]     bc_def   `BCDefinition` from which operators are created and added
202*77aa5ad7SJames Wright   @param[in,out] op_ifunc Composite operator for IFunction operators to be added to
203*77aa5ad7SJames Wright   @param[in,out] op_ijac  Composite operator for IJacobian operators to be added to
204*77aa5ad7SJames Wright **/
205*77aa5ad7SJames Wright PetscErrorCode BCDefinitionAddOperators(BCDefinition bc_def, CeedOperator op_ifunc, CeedOperator op_ijac) {
206*77aa5ad7SJames Wright   Ceed            ceed = CeedOperatorReturnCeed(op_ifunc);
207*77aa5ad7SJames Wright   CeedQFunction   qf_ifunction, qf_ijacobian;
208*77aa5ad7SJames Wright   DMLabel         face_sets_label;
209*77aa5ad7SJames Wright   PetscInt        num_face_set_values;
210*77aa5ad7SJames Wright   const PetscInt *face_set_values;
211*77aa5ad7SJames Wright 
212*77aa5ad7SJames Wright   PetscFunctionBeginUser;
213*77aa5ad7SJames Wright   if (!bc_def->CreateIFunctionQF || !bc_def->AddIFunctionOperator) PetscFunctionReturn(PETSC_SUCCESS);
214*77aa5ad7SJames Wright   PetscBool add_ijac = (!bc_def->CreateIJacobianQF || !bc_def->AddIJacobianOperator) ? PETSC_FALSE : PETSC_TRUE;
215*77aa5ad7SJames Wright   PetscCheck(bc_def->dm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "BCDefinition must have DM set using BCDefinitionSetDM()");
216*77aa5ad7SJames Wright 
217*77aa5ad7SJames Wright   PetscCall(bc_def->CreateIFunctionQF(bc_def, &qf_ifunction));
218*77aa5ad7SJames Wright   if (add_ijac) PetscCall(bc_def->CreateIJacobianQF(bc_def, &qf_ijacobian));
219*77aa5ad7SJames Wright 
220*77aa5ad7SJames Wright   PetscCall(DMGetLabel(bc_def->dm, "Face Sets", &face_sets_label));
221*77aa5ad7SJames Wright   PetscCall(BCDefinitionGetInfo(bc_def, NULL, &num_face_set_values, &face_set_values));
222*77aa5ad7SJames Wright   for (PetscInt f = 0; f < num_face_set_values; f++) {
223*77aa5ad7SJames Wright     DMLabel  face_orientation_label;
224*77aa5ad7SJames Wright     PetscInt num_orientations_values, *orientation_values;
225*77aa5ad7SJames Wright 
226*77aa5ad7SJames Wright     {
227*77aa5ad7SJames Wright       char *face_orientation_label_name;
228*77aa5ad7SJames Wright 
229*77aa5ad7SJames Wright       PetscCall(DMPlexCreateFaceLabel(bc_def->dm, face_set_values[f], &face_orientation_label_name));
230*77aa5ad7SJames Wright       PetscCall(DMGetLabel(bc_def->dm, face_orientation_label_name, &face_orientation_label));
231*77aa5ad7SJames Wright       PetscCall(PetscFree(face_orientation_label_name));
232*77aa5ad7SJames Wright     }
233*77aa5ad7SJames Wright     PetscCall(DMLabelCreateGlobalValueArray(bc_def->dm, face_orientation_label, &num_orientations_values, &orientation_values));
234*77aa5ad7SJames Wright     for (PetscInt o = 0; o < num_orientations_values; o++) {
235*77aa5ad7SJames Wright       CeedOperator sub_op_ifunc;
236*77aa5ad7SJames Wright       PetscInt     orientation = orientation_values[o];
237*77aa5ad7SJames Wright 
238*77aa5ad7SJames Wright       PetscCall(bc_def->AddIFunctionOperator(bc_def, face_orientation_label, orientation, qf_ifunction, op_ifunc, &sub_op_ifunc));
239*77aa5ad7SJames Wright       if (add_ijac) PetscCall(bc_def->AddIJacobianOperator(bc_def, sub_op_ifunc, face_orientation_label, orientation, qf_ijacobian, op_ijac));
240*77aa5ad7SJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&sub_op_ifunc));
241*77aa5ad7SJames Wright     }
242*77aa5ad7SJames Wright     PetscCall(PetscFree(orientation_values));
243*77aa5ad7SJames Wright   }
244*77aa5ad7SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction));
245*77aa5ad7SJames Wright   if (add_ijac) PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian));
246*77aa5ad7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
247*77aa5ad7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
248*77aa5ad7SJames Wright }
249