1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 /// @file 4 /// Functions for setting up and projecting the velocity gradient 5 6 #include "../qfunctions/velocity_gradient_projection.h" 7 8 #include <petscdmplex.h> 9 10 #include <navierstokes.h> 11 12 PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) { 13 PetscSection section; 14 15 PetscFunctionBeginUser; 16 grad_velo_proj->num_comp = 9; // 9 velocity gradient 17 18 PetscCall(DMClone(user->dm, &grad_velo_proj->dm)); 19 PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection")); 20 21 PetscCall( 22 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)); 23 24 PetscCall(DMGetLocalSection(grad_velo_proj->dm, §ion)); 25 PetscCall(PetscSectionSetFieldName(section, 0, "")); 26 PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX")); 27 PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY")); 28 PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ")); 29 PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX")); 30 PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY")); 31 PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ")); 32 PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX")); 33 PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY")); 34 PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ")); 35 PetscFunctionReturn(PETSC_SUCCESS); 36 }; 37 38 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input, 39 CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj) { 40 NodalProjectionData grad_velo_proj; 41 CeedOperator op_rhs_assemble, op_mass; 42 CeedQFunction qf_rhs_assemble, qf_mass; 43 CeedBasis basis_grad_velo; 44 CeedElemRestriction elem_restr_grad_velo; 45 PetscInt dim, label_value = 0, height = 0, dm_field = 0; 46 CeedInt num_comp_x, num_comp_input, q_data_size; 47 DMLabel domain_label = NULL; 48 49 PetscFunctionBeginUser; 50 PetscCall(PetscNew(&grad_velo_proj)); 51 52 PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree)); 53 54 // -- Get Pre-requisite things 55 PetscCall(DMGetDimension(grad_velo_proj->dm, &dim)); 56 PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x)); 57 PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_input)); 58 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size)); 59 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo)); 60 61 PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo)); 62 63 // -- Build RHS operator 64 switch (state_var_input) { 65 case STATEVAR_PRIMITIVE: 66 PetscCallCeed( 67 ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble)); 68 break; 69 case STATEVAR_CONSERVATIVE: 70 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, 71 &qf_rhs_assemble)); 72 break; 73 case STATEVAR_ENTROPY: 74 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Entropy, VelocityGradientProjectionRHS_Entropy_loc, 75 &qf_rhs_assemble)); 76 break; 77 } 78 79 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context)); 80 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP)); 81 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD)); 82 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE)); 83 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP)); 84 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP)); 85 86 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble)); 87 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 88 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 89 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 90 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 91 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 92 93 PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx)); 94 95 // -- Build Mass operator 96 PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass)); 97 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 98 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 99 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 100 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 101 102 { // -- Setup KSP for L^2 projection with lumped mass operator 103 Mat mat_mass; 104 MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm); 105 106 PetscCall(MatCreateCeed(grad_velo_proj->dm, grad_velo_proj->dm, op_mass, NULL, &mat_mass)); 107 108 PetscCall(KSPCreate(comm, &grad_velo_proj->ksp)); 109 PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_")); 110 { 111 PC pc; 112 PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc)); 113 PetscCall(PCSetType(pc, PCJACOBI)); 114 PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 115 PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY)); 116 } 117 PetscCall(KSPSetFromOptions_WithMatCeed(grad_velo_proj->ksp, mat_mass)); 118 PetscCall(MatDestroy(&mat_mass)); 119 } 120 121 *pgrad_velo_proj = grad_velo_proj; 122 123 PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo)); 124 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 125 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble)); 126 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 127 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble)); 128 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 129 PetscFunctionReturn(PETSC_SUCCESS); 130 } 131 132 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) { 133 OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx; 134 135 PetscFunctionBeginUser; 136 PetscCall(PetscLogEventBegin(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 137 PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx)); 138 139 PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient)); 140 PetscCall(PetscLogEventEnd(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 141 PetscFunctionReturn(PETSC_SUCCESS); 142 } 143