15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2999ff5c7SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3999ff5c7SJames Wright // 4999ff5c7SJames Wright // SPDX-License-Identifier: BSD-2-Clause 5999ff5c7SJames Wright // 6999ff5c7SJames Wright // This file is part of CEED: http://github.com/ceed 7999ff5c7SJames Wright /// @file 8999ff5c7SJames Wright /// Functions for setting up and projecting the velocity gradient 9999ff5c7SJames Wright 10999ff5c7SJames Wright #include "../qfunctions/velocity_gradient_projection.h" 11999ff5c7SJames Wright 12999ff5c7SJames Wright #include <petscdmplex.h> 13999ff5c7SJames Wright 14999ff5c7SJames Wright #include "../navierstokes.h" 15999ff5c7SJames Wright 16999ff5c7SJames Wright PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) { 17999ff5c7SJames Wright PetscSection section; 18999ff5c7SJames Wright 19999ff5c7SJames Wright PetscFunctionBeginUser; 20999ff5c7SJames Wright grad_velo_proj->num_comp = 9; // 9 velocity gradient 21999ff5c7SJames Wright 22999ff5c7SJames Wright PetscCall(DMClone(user->dm, &grad_velo_proj->dm)); 23999ff5c7SJames Wright PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection")); 24999ff5c7SJames Wright 25d68a66c4SJames Wright PetscCall( 26d68a66c4SJames 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)); 27999ff5c7SJames Wright 28999ff5c7SJames Wright PetscCall(DMGetLocalSection(grad_velo_proj->dm, §ion)); 29999ff5c7SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 30999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX")); 31999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY")); 32999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ")); 33999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX")); 34999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY")); 35999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ")); 36999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX")); 37999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY")); 38999ff5c7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ")); 39ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 40999ff5c7SJames Wright }; 413909d146SJames Wright 42*731c13d7SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input, 437618d7b7SJames Wright CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj) { 447618d7b7SJames Wright NodalProjectionData grad_velo_proj; 453909d146SJames Wright CeedOperator op_rhs_assemble, op_mass; 463909d146SJames Wright CeedQFunction qf_rhs_assemble, qf_mass; 473909d146SJames Wright CeedBasis basis_grad_velo; 483909d146SJames Wright CeedElemRestriction elem_restr_grad_velo; 49bb85d312SJames Wright PetscInt dim, label_value = 0, height = 0, dm_field = 0; 507618d7b7SJames Wright CeedInt num_comp_x, num_comp_input, q_data_size; 51bb85d312SJames Wright DMLabel domain_label = NULL; 523909d146SJames Wright 533909d146SJames Wright PetscFunctionBeginUser; 547618d7b7SJames Wright PetscCall(PetscNew(&grad_velo_proj)); 553909d146SJames Wright 563909d146SJames Wright PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree)); 573909d146SJames Wright 583909d146SJames Wright // -- Get Pre-requisite things 593909d146SJames Wright PetscCall(DMGetDimension(grad_velo_proj->dm, &dim)); 60a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x)); 617618d7b7SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_input)); 62a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size)); 63bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo)); 643909d146SJames Wright 65bb85d312SJames Wright PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo)); 663909d146SJames Wright 673909d146SJames Wright // -- Build RHS operator 687618d7b7SJames Wright switch (state_var_input) { 693909d146SJames Wright case STATEVAR_PRIMITIVE: 70a424bcd0SJames Wright PetscCallCeed( 71a424bcd0SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble)); 723909d146SJames Wright break; 733909d146SJames Wright case STATEVAR_CONSERVATIVE: 74a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, 75a424bcd0SJames Wright &qf_rhs_assemble)); 763909d146SJames Wright break; 773909d146SJames Wright default: 783909d146SJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable"); 793909d146SJames Wright } 803909d146SJames Wright 81a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context)); 827618d7b7SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP)); 837618d7b7SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD)); 84a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE)); 85a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP)); 86a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP)); 873909d146SJames Wright 88a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble)); 897618d7b7SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 907618d7b7SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE)); 91356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 92a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 93a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 943909d146SJames Wright 95e66242d1SJames Wright PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx)); 963909d146SJames Wright 973909d146SJames Wright // -- Build Mass operator 983909d146SJames Wright PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass)); 99a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 100a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 101356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 102a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE)); 1033909d146SJames Wright 1043909d146SJames Wright { // -- Setup KSP for L^2 projection with lumped mass operator 1057f2a9303SJames Wright Mat mat_mass; 1063909d146SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm); 1073909d146SJames Wright 1087f2a9303SJames Wright PetscCall(MatCeedCreate(grad_velo_proj->dm, grad_velo_proj->dm, op_mass, NULL, &mat_mass)); 1093909d146SJames Wright 1103909d146SJames Wright PetscCall(KSPCreate(comm, &grad_velo_proj->ksp)); 1113909d146SJames Wright PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_")); 1123909d146SJames Wright { 1133909d146SJames Wright PC pc; 1143909d146SJames Wright PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc)); 1153909d146SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 1163909d146SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 1173909d146SJames Wright PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY)); 1183909d146SJames Wright } 1197f2a9303SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(grad_velo_proj->ksp, mat_mass)); 120cde30410SJames Wright PetscCall(MatDestroy(&mat_mass)); 1213909d146SJames Wright } 1223909d146SJames Wright 1237618d7b7SJames Wright *pgrad_velo_proj = grad_velo_proj; 1247618d7b7SJames Wright 125a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo)); 126a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 127a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble)); 128a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 129a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble)); 130a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 131ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1323909d146SJames Wright } 133f84dcd56SJames Wright 1347618d7b7SJames Wright PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) { 135f84dcd56SJames Wright OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx; 136f84dcd56SJames Wright 137f84dcd56SJames Wright PetscFunctionBeginUser; 1383737f832SJames Wright PetscCall(PetscLogEventBegin(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 139f84dcd56SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx)); 140f84dcd56SJames Wright 141f84dcd56SJames Wright PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient)); 1423737f832SJames Wright PetscCall(PetscLogEventEnd(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0)); 143ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 144f84dcd56SJames Wright } 145