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 120c373b74SJames Wright 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 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; 4515c18037SJames Wright PetscInt dim, label_value = 0, height = 0, dm_field = 0; 46f6ac214eSJames Wright CeedInt num_comp_x, num_comp_input, q_data_size; 4715c18037SJames Wright DMLabel domain_label = NULL; 489ab09d51SJames Wright 499ab09d51SJames Wright PetscFunctionBeginUser; 50f6ac214eSJames Wright PetscCall(PetscNew(&grad_velo_proj)); 510c373b74SJames Wright PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, honee, honee->app_ctx->degree)); 529ab09d51SJames Wright 539ab09d51SJames Wright // -- Get Pre-requisite things 549ab09d51SJames Wright PetscCall(DMGetDimension(grad_velo_proj->dm, &dim)); 55e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_x, &num_comp_x)); 56e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_input)); 5715c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo)); 5815c18037SJames Wright PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo)); 59e3663b90SJames Wright PetscCall(QDataGet(ceed, grad_velo_proj->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, 60e3663b90SJames Wright &q_data, &q_data_size)); 619ab09d51SJames Wright 624aea4664SJames Wright { // -- Build RHS operator 634aea4664SJames Wright CeedOperator op_rhs_assemble; 644aea4664SJames Wright CeedQFunction qf_rhs_assemble; 654aea4664SJames Wright 66f6ac214eSJames Wright switch (state_var_input) { 679ab09d51SJames Wright case STATEVAR_PRIMITIVE: 68b4c37c5cSJames Wright PetscCallCeed( 69b4c37c5cSJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble)); 709ab09d51SJames Wright break; 719ab09d51SJames Wright case STATEVAR_CONSERVATIVE: 72b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, 73b4c37c5cSJames Wright &qf_rhs_assemble)); 749ab09d51SJames Wright break; 759b103f75SJames Wright case STATEVAR_ENTROPY: 769b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Entropy, VelocityGradientProjectionRHS_Entropy_loc, 779b103f75SJames Wright &qf_rhs_assemble)); 789b103f75SJames Wright break; 799ab09d51SJames Wright } 809ab09d51SJames Wright 81e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfctx)); 82f6ac214eSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP)); 83f6ac214eSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD)); 84b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE)); 85b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP)); 86b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP)); 879ab09d51SJames Wright 88b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble)); 89f6ac214eSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 90f6ac214eSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 91be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 92e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 93b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 949ab09d51SJames Wright 950c373b74SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx)); 969ab09d51SJames Wright 974aea4664SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble)); 984aea4664SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble)); 994aea4664SJames Wright } 1004aea4664SJames Wright 1014aea4664SJames Wright { // -- Build Mass operator 1024aea4664SJames Wright CeedOperator op_mass; 1034aea4664SJames Wright CeedQFunction qf_mass; 1044aea4664SJames Wright Mat mat_mass; 1054aea4664SJames Wright 106*64dd23feSJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass)); 107b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 108b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 109be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 110b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 1119ab09d51SJames Wright 112000d2032SJeremy L Thompson PetscCall(MatCreateCeed(grad_velo_proj->dm, grad_velo_proj->dm, op_mass, NULL, &mat_mass)); 1139ab09d51SJames Wright 1144aea4664SJames Wright // Setup KSP for L^2 projection with lumped mass operator 1154aea4664SJames Wright PetscCall(KSPCreate(PetscObjectComm((PetscObject)grad_velo_proj->dm), &grad_velo_proj->ksp)); 1169ab09d51SJames Wright PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_")); 1179ab09d51SJames Wright { 1189ab09d51SJames Wright PC pc; 1199ab09d51SJames Wright PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc)); 1209ab09d51SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 1219ab09d51SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 1229ab09d51SJames Wright PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY)); 1239ab09d51SJames Wright } 1243170c09fSJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(grad_velo_proj->ksp, mat_mass)); 1254aea4664SJames Wright 1264aea4664SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 1274aea4664SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 128b8462d76SJames Wright PetscCall(MatDestroy(&mat_mass)); 1299ab09d51SJames Wright } 1309ab09d51SJames Wright 131f6ac214eSJames Wright *pgrad_velo_proj = grad_velo_proj; 132f6ac214eSJames Wright 133b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo)); 134be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 135be29160dSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 136b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 137d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1389ab09d51SJames Wright } 139a2401be5SJames Wright 140f6ac214eSJames Wright PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) { 141a2401be5SJames Wright OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx; 142a2401be5SJames Wright 143a2401be5SJames Wright PetscFunctionBeginUser; 144855536edSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 145a2401be5SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx)); 146a2401be5SJames Wright 147a2401be5SJames Wright PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient)); 148855536edSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 149d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 150a2401be5SJames Wright } 151