1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3457a5831SJames Wright /// @file
4457a5831SJames Wright /// Functions for setting up and projecting the velocity gradient
5457a5831SJames Wright
6457a5831SJames Wright #include "../qfunctions/velocity_gradient_projection.h"
7457a5831SJames Wright
8457a5831SJames Wright #include <petscdmplex.h>
9457a5831SJames Wright
10149fb536SJames Wright #include <navierstokes.h>
11457a5831SJames Wright
VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj,Honee honee,PetscInt degree)12cf8f12c8SJames Wright static PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, Honee honee, PetscInt degree) {
13457a5831SJames Wright PetscSection section;
14457a5831SJames Wright
15457a5831SJames Wright PetscFunctionBeginUser;
16457a5831SJames Wright grad_velo_proj->num_comp = 9; // 9 velocity gradient
17457a5831SJames Wright
180c373b74SJames Wright PetscCall(DMClone(honee->dm, &grad_velo_proj->dm));
190dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(grad_velo_proj->dm, PETSC_TRUE));
20457a5831SJames Wright PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
21457a5831SJames Wright
220c373b74SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, 1, &grad_velo_proj->num_comp,
230c373b74SJames Wright grad_velo_proj->dm));
24457a5831SJames Wright
25457a5831SJames Wright PetscCall(DMGetLocalSection(grad_velo_proj->dm, §ion));
26457a5831SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, ""));
27457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
28457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
29457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
30457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
31457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
32457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
33457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
34457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
35457a5831SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
36d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
37457a5831SJames Wright };
389ab09d51SJames Wright
VelocityGradientProjectionSetup(Ceed ceed,Honee honee,ProblemData problem,StateVariable state_var_input,CeedElemRestriction elem_restr_input,CeedBasis basis_input,NodalProjectionData * pgrad_velo_proj)39e3663b90SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input,
40f6ac214eSJames Wright CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj) {
41f6ac214eSJames Wright NodalProjectionData grad_velo_proj;
429ab09d51SJames Wright CeedBasis basis_grad_velo;
43be29160dSJames Wright CeedElemRestriction elem_restr_grad_velo, elem_restr_qd;
44be29160dSJames Wright CeedVector q_data;
45*4fe35dceSJames Wright PetscInt dim, height = 0, dm_field = 0, num_comp_input;
46cf8f12c8SJames Wright CeedInt q_data_size;
479ab09d51SJames Wright
489ab09d51SJames Wright PetscFunctionBeginUser;
49f6ac214eSJames Wright PetscCall(PetscNew(&grad_velo_proj));
500c373b74SJames Wright PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, honee, honee->app_ctx->degree));
519ab09d51SJames Wright
529ab09d51SJames Wright // -- Get Pre-requisite things
539ab09d51SJames Wright PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
54cf8f12c8SJames Wright PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_input));
55e3db12f8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field,
56e3db12f8SJames Wright &elem_restr_grad_velo));
57e3db12f8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, grad_velo_proj->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_grad_velo));
589018c49aSJames Wright PetscCall(QDataGet(ceed, grad_velo_proj->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &q_data, &q_data_size));
599ab09d51SJames Wright
604aea4664SJames Wright { // -- Build RHS operator
614aea4664SJames Wright CeedOperator op_rhs_assemble;
624aea4664SJames Wright CeedQFunction qf_rhs_assemble;
63*4fe35dceSJames Wright CeedElemRestriction elem_restr_x;
64*4fe35dceSJames Wright CeedBasis basis_x;
65*4fe35dceSJames Wright CeedVector x_coord;
66*4fe35dceSJames Wright PetscInt num_comp_x;
674aea4664SJames Wright
68*4fe35dceSJames Wright PetscCall(DMGetCoordinateNumComps(honee->dm, &num_comp_x));
69*4fe35dceSJames Wright PetscCall(DMPlexCeedCoordinateCreateField(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, &elem_restr_x, &basis_x, &x_coord));
70f6ac214eSJames Wright switch (state_var_input) {
719ab09d51SJames Wright case STATEVAR_PRIMITIVE:
729eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc,
739eadbee4SJames Wright &qf_rhs_assemble));
749ab09d51SJames Wright break;
759ab09d51SJames Wright case STATEVAR_CONSERVATIVE:
76b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc,
77b4c37c5cSJames Wright &qf_rhs_assemble));
789ab09d51SJames Wright break;
799b103f75SJames Wright case STATEVAR_ENTROPY:
809b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Entropy, VelocityGradientProjectionRHS_Entropy_loc,
819b103f75SJames Wright &qf_rhs_assemble));
829b103f75SJames Wright break;
839ab09d51SJames Wright }
849ab09d51SJames Wright
85e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfctx));
86f6ac214eSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP));
87f6ac214eSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD));
88b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE));
89b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP));
90b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP));
919ab09d51SJames Wright
92b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble));
93f6ac214eSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
94f6ac214eSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
95be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
96*4fe35dceSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", elem_restr_x, basis_x, x_coord));
97b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
989ab09d51SJames Wright
990c373b74SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
1009ab09d51SJames Wright
101*4fe35dceSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x));
102*4fe35dceSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x));
103*4fe35dceSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&x_coord));
1044aea4664SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble));
1054aea4664SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble));
1064aea4664SJames Wright }
1074aea4664SJames Wright
1084aea4664SJames Wright { // -- Build Mass operator
1094aea4664SJames Wright CeedOperator op_mass;
1104aea4664SJames Wright CeedQFunction qf_mass;
1114aea4664SJames Wright Mat mat_mass;
1124aea4664SJames Wright
11364dd23feSJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
114b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
115b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
116be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
117b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
1189ab09d51SJames Wright
119000d2032SJeremy L Thompson PetscCall(MatCreateCeed(grad_velo_proj->dm, grad_velo_proj->dm, op_mass, NULL, &mat_mass));
1209ab09d51SJames Wright
1214aea4664SJames Wright // Setup KSP for L^2 projection with lumped mass operator
1224aea4664SJames Wright PetscCall(KSPCreate(PetscObjectComm((PetscObject)grad_velo_proj->dm), &grad_velo_proj->ksp));
1239ab09d51SJames Wright PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
1249ab09d51SJames Wright {
1259ab09d51SJames Wright PC pc;
1269ab09d51SJames Wright PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
1279ab09d51SJames Wright PetscCall(PCSetType(pc, PCJACOBI));
1289ab09d51SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1299ab09d51SJames Wright PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
1309ab09d51SJames Wright }
1313170c09fSJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(grad_velo_proj->ksp, mat_mass));
1324aea4664SJames Wright
1334aea4664SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
1344aea4664SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
135b8462d76SJames Wright PetscCall(MatDestroy(&mat_mass));
1369ab09d51SJames Wright }
1379ab09d51SJames Wright
138f6ac214eSJames Wright *pgrad_velo_proj = grad_velo_proj;
139f6ac214eSJames Wright
140b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo));
141be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
142be29160dSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
143b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
144d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
1459ab09d51SJames Wright }
146a2401be5SJames Wright
VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj,Vec Q_loc,Vec VelocityGradient)147f6ac214eSJames Wright PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) {
148a2401be5SJames Wright OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx;
149a2401be5SJames Wright
150a2401be5SJames Wright PetscFunctionBeginUser;
151ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
152a2401be5SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
153a2401be5SJames Wright
154a2401be5SJames Wright PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
155ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
156d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
157a2401be5SJames Wright }
158