xref: /honee/src/velocity_gradient_projection.c (revision ae2b091fac884a554e48acc4b4c187524c2a2818)
1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2*ae2b091fSJames 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));
19457a5831SJames Wright   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
20457a5831SJames Wright 
21da4ca0cfSJames Wright   PetscCall(
22da4ca0cfSJames 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));
23457a5831SJames Wright 
24457a5831SJames Wright   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
25457a5831SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
26457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
27457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
28457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
29457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
30457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
31457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
32457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
33457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
34457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
35d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
36457a5831SJames Wright };
379ab09d51SJames Wright 
38991aef52SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input,
39f6ac214eSJames Wright                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj) {
40f6ac214eSJames Wright   NodalProjectionData grad_velo_proj;
419ab09d51SJames Wright   CeedOperator        op_rhs_assemble, op_mass;
429ab09d51SJames Wright   CeedQFunction       qf_rhs_assemble, qf_mass;
439ab09d51SJames Wright   CeedBasis           basis_grad_velo;
449ab09d51SJames Wright   CeedElemRestriction elem_restr_grad_velo;
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));
519ab09d51SJames Wright 
529ab09d51SJames Wright   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
539ab09d51SJames Wright 
549ab09d51SJames Wright   // -- Get Pre-requisite things
559ab09d51SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
56b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
57f6ac214eSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_input));
58b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
5915c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo));
609ab09d51SJames Wright 
6115c18037SJames Wright   PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo));
629ab09d51SJames Wright 
639ab09d51SJames Wright   // -- Build RHS operator
64f6ac214eSJames Wright   switch (state_var_input) {
659ab09d51SJames Wright     case STATEVAR_PRIMITIVE:
66b4c37c5cSJames Wright       PetscCallCeed(
67b4c37c5cSJames Wright           ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble));
689ab09d51SJames Wright       break;
699ab09d51SJames Wright     case STATEVAR_CONSERVATIVE:
70b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc,
71b4c37c5cSJames Wright                                                       &qf_rhs_assemble));
729ab09d51SJames Wright       break;
739b103f75SJames Wright     case STATEVAR_ENTROPY:
749b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Entropy, VelocityGradientProjectionRHS_Entropy_loc,
759b103f75SJames Wright                                                       &qf_rhs_assemble));
769b103f75SJames Wright       break;
779ab09d51SJames Wright   }
789ab09d51SJames Wright 
79b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context));
80f6ac214eSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP));
81f6ac214eSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD));
82b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE));
83b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP));
84b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP));
859ab09d51SJames Wright 
86b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble));
87f6ac214eSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
88f6ac214eSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
8958e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
90b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
91b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
929ab09d51SJames Wright 
9349f5db74SJames Wright   PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
949ab09d51SJames Wright 
959ab09d51SJames Wright   // -- Build Mass operator
969ab09d51SJames Wright   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
97b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
98b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
9958e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
100b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
1019ab09d51SJames Wright 
1029ab09d51SJames Wright   {  // -- Setup KSP for L^2 projection with lumped mass operator
1033170c09fSJames Wright     Mat      mat_mass;
1049ab09d51SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
1059ab09d51SJames Wright 
106000d2032SJeremy L Thompson     PetscCall(MatCreateCeed(grad_velo_proj->dm, grad_velo_proj->dm, op_mass, NULL, &mat_mass));
1079ab09d51SJames Wright 
1089ab09d51SJames Wright     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
1099ab09d51SJames Wright     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
1109ab09d51SJames Wright     {
1119ab09d51SJames Wright       PC pc;
1129ab09d51SJames Wright       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
1139ab09d51SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
1149ab09d51SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1159ab09d51SJames Wright       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
1169ab09d51SJames Wright     }
1173170c09fSJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(grad_velo_proj->ksp, mat_mass));
118b8462d76SJames Wright     PetscCall(MatDestroy(&mat_mass));
1199ab09d51SJames Wright   }
1209ab09d51SJames Wright 
121f6ac214eSJames Wright   *pgrad_velo_proj = grad_velo_proj;
122f6ac214eSJames Wright 
123b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo));
124b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
125b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble));
126b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
127b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble));
128b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
129d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1309ab09d51SJames Wright }
131a2401be5SJames Wright 
132f6ac214eSJames Wright PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) {
133a2401be5SJames Wright   OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx;
134a2401be5SJames Wright 
135a2401be5SJames Wright   PetscFunctionBeginUser;
136855536edSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
137a2401be5SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
138a2401be5SJames Wright 
139a2401be5SJames Wright   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
140855536edSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
141d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
142a2401be5SJames Wright }
143