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 12457a5831SJames Wright PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) { 13457a5831SJames Wright PetscSection section; 14457a5831SJames Wright 15457a5831SJames Wright PetscFunctionBeginUser; 16457a5831SJames Wright grad_velo_proj->num_comp = 9; // 9 velocity gradient 17457a5831SJames Wright 18457a5831SJames Wright PetscCall(DMClone(user->dm, &grad_velo_proj->dm)); 19*0dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(grad_velo_proj->dm, PETSC_TRUE)); 20457a5831SJames Wright PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection")); 21457a5831SJames Wright 22da4ca0cfSJames Wright PetscCall( 23da4ca0cfSJames Wright DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &grad_velo_proj->num_comp, 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 39991aef52SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, 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 CeedOperator op_rhs_assemble, op_mass; 439ab09d51SJames Wright CeedQFunction qf_rhs_assemble, qf_mass; 449ab09d51SJames Wright CeedBasis basis_grad_velo; 459ab09d51SJames Wright CeedElemRestriction elem_restr_grad_velo; 4615c18037SJames Wright PetscInt dim, label_value = 0, height = 0, dm_field = 0; 47f6ac214eSJames Wright CeedInt num_comp_x, num_comp_input, q_data_size; 4815c18037SJames Wright DMLabel domain_label = NULL; 499ab09d51SJames Wright 509ab09d51SJames Wright PetscFunctionBeginUser; 51f6ac214eSJames Wright PetscCall(PetscNew(&grad_velo_proj)); 529ab09d51SJames Wright 539ab09d51SJames Wright PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree)); 549ab09d51SJames Wright 559ab09d51SJames Wright // -- Get Pre-requisite things 569ab09d51SJames Wright PetscCall(DMGetDimension(grad_velo_proj->dm, &dim)); 57b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x)); 58f6ac214eSJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_input)); 59b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size)); 6015c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo)); 619ab09d51SJames Wright 6215c18037SJames Wright PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo)); 639ab09d51SJames Wright 649ab09d51SJames Wright // -- Build RHS operator 65f6ac214eSJames Wright switch (state_var_input) { 669ab09d51SJames Wright case STATEVAR_PRIMITIVE: 67b4c37c5cSJames Wright PetscCallCeed( 68b4c37c5cSJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble)); 699ab09d51SJames Wright break; 709ab09d51SJames Wright case STATEVAR_CONSERVATIVE: 71b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, 72b4c37c5cSJames Wright &qf_rhs_assemble)); 739ab09d51SJames Wright break; 749b103f75SJames Wright case STATEVAR_ENTROPY: 759b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Entropy, VelocityGradientProjectionRHS_Entropy_loc, 769b103f75SJames Wright &qf_rhs_assemble)); 779b103f75SJames Wright break; 789ab09d51SJames Wright } 799ab09d51SJames Wright 80e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfctx)); 81f6ac214eSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP)); 82f6ac214eSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD)); 83b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE)); 84b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP)); 85b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP)); 869ab09d51SJames Wright 87b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble)); 88f6ac214eSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 89f6ac214eSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 9058e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 91b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 92b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 939ab09d51SJames Wright 9449f5db74SJames Wright PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx)); 959ab09d51SJames Wright 969ab09d51SJames Wright // -- Build Mass operator 979ab09d51SJames Wright PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass)); 98b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 99b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 10058e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 101b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 1029ab09d51SJames Wright 1039ab09d51SJames Wright { // -- Setup KSP for L^2 projection with lumped mass operator 1043170c09fSJames Wright Mat mat_mass; 1059ab09d51SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm); 1069ab09d51SJames Wright 107000d2032SJeremy L Thompson PetscCall(MatCreateCeed(grad_velo_proj->dm, grad_velo_proj->dm, op_mass, NULL, &mat_mass)); 1089ab09d51SJames Wright 1099ab09d51SJames Wright PetscCall(KSPCreate(comm, &grad_velo_proj->ksp)); 1109ab09d51SJames Wright PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_")); 1119ab09d51SJames Wright { 1129ab09d51SJames Wright PC pc; 1139ab09d51SJames Wright PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc)); 1149ab09d51SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 1159ab09d51SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 1169ab09d51SJames Wright PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY)); 1179ab09d51SJames Wright } 1183170c09fSJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(grad_velo_proj->ksp, mat_mass)); 119b8462d76SJames Wright PetscCall(MatDestroy(&mat_mass)); 1209ab09d51SJames Wright } 1219ab09d51SJames Wright 122f6ac214eSJames Wright *pgrad_velo_proj = grad_velo_proj; 123f6ac214eSJames Wright 124b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo)); 125b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 126b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble)); 127b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 128b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble)); 129b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 130d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1319ab09d51SJames Wright } 132a2401be5SJames Wright 133f6ac214eSJames Wright PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) { 134a2401be5SJames Wright OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx; 135a2401be5SJames Wright 136a2401be5SJames Wright PetscFunctionBeginUser; 137855536edSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 138a2401be5SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx)); 139a2401be5SJames Wright 140a2401be5SJames Wright PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient)); 141855536edSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 142d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 143a2401be5SJames Wright } 144