xref: /honee/src/velocity_gradient_projection.c (revision a2401be548428d9c79d7c84f79af63730cee72ad)
1457a5831SJames Wright // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
2457a5831SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3457a5831SJames Wright //
4457a5831SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5457a5831SJames Wright //
6457a5831SJames Wright // This file is part of CEED:  http://github.com/ceed
7457a5831SJames Wright /// @file
8457a5831SJames Wright /// Functions for setting up and projecting the velocity gradient
9457a5831SJames Wright 
10457a5831SJames Wright #include "../qfunctions/velocity_gradient_projection.h"
11457a5831SJames Wright 
12457a5831SJames Wright #include <petscdmplex.h>
13457a5831SJames Wright 
14457a5831SJames Wright #include "../navierstokes.h"
15457a5831SJames Wright 
16457a5831SJames Wright PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) {
17457a5831SJames Wright   PetscFE      fe;
18457a5831SJames Wright   PetscSection section;
19457a5831SJames Wright   PetscInt     dim;
20457a5831SJames Wright 
21457a5831SJames Wright   PetscFunctionBeginUser;
22457a5831SJames Wright   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
23457a5831SJames Wright 
24457a5831SJames Wright   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
25457a5831SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
26457a5831SJames Wright   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
27457a5831SJames Wright 
28457a5831SJames Wright   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, grad_velo_proj->num_comp, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
29457a5831SJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "Velocity Gradient Projection"));
30457a5831SJames Wright   PetscCall(DMAddField(grad_velo_proj->dm, NULL, (PetscObject)fe));
31457a5831SJames Wright   PetscCall(DMCreateDS(grad_velo_proj->dm));
32457a5831SJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(grad_velo_proj->dm, PETSC_DETERMINE, NULL));
33457a5831SJames Wright 
34457a5831SJames Wright   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
35457a5831SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
36457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
37457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
38457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
39457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
40457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
41457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
42457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
43457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
44457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
45457a5831SJames Wright 
46457a5831SJames Wright   PetscCall(PetscFEDestroy(&fe));
47457a5831SJames Wright   PetscFunctionReturn(0);
48457a5831SJames Wright };
499ab09d51SJames Wright 
509ab09d51SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
519ab09d51SJames Wright   OperatorApplyContext mass_matop_ctx;
529ab09d51SJames Wright   CeedOperator         op_rhs_assemble, op_mass;
539ab09d51SJames Wright   CeedQFunction        qf_rhs_assemble, qf_mass;
549ab09d51SJames Wright   CeedBasis            basis_grad_velo;
559ab09d51SJames Wright   CeedVector           q_ceed, rhs_ceed, mass_output;
569ab09d51SJames Wright   CeedElemRestriction  elem_restr_grad_velo;
579ab09d51SJames Wright   PetscInt             dim, num_comp_x, num_comp_q, q_data_size, num_qpts_1d, num_nodes_1d;
589ab09d51SJames Wright 
599ab09d51SJames Wright   PetscFunctionBeginUser;
609ab09d51SJames Wright   PetscCall(PetscNew(&user->grad_velo_proj));
619ab09d51SJames Wright   NodalProjectionData grad_velo_proj = user->grad_velo_proj;
629ab09d51SJames Wright 
639ab09d51SJames Wright   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
649ab09d51SJames Wright 
659ab09d51SJames Wright   // -- Get Pre-requisite things
669ab09d51SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
679ab09d51SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
689ab09d51SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
699ab09d51SJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &num_qpts_1d);
709ab09d51SJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &num_nodes_1d);
719ab09d51SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
729ab09d51SJames Wright 
739ab09d51SJames Wright   PetscCall(GetRestrictionForDomain(ceed, grad_velo_proj->dm, 0, 0, 0, num_qpts_1d, q_data_size, &elem_restr_grad_velo, NULL, NULL));
749ab09d51SJames Wright 
759ab09d51SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, grad_velo_proj->num_comp, num_nodes_1d, num_qpts_1d, CEED_GAUSS, &basis_grad_velo);
769ab09d51SJames Wright 
779ab09d51SJames Wright   // -- Build RHS operator
789ab09d51SJames Wright   switch (user->phys->state_var) {
799ab09d51SJames Wright     case STATEVAR_PRIMITIVE:
809ab09d51SJames Wright       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble);
819ab09d51SJames Wright       break;
829ab09d51SJames Wright     case STATEVAR_CONSERVATIVE:
839ab09d51SJames Wright       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, &qf_rhs_assemble);
849ab09d51SJames Wright       break;
859ab09d51SJames Wright     default:
869ab09d51SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
879ab09d51SJames Wright   }
889ab09d51SJames Wright 
899ab09d51SJames Wright   CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context);
909ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_q, CEED_EVAL_INTERP);
919ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
929ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE);
939ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP);
949ab09d51SJames Wright   CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP);
959ab09d51SJames Wright 
969ab09d51SJames Wright   CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble);
979ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
989ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
999ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
1009ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
1019ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
1029ab09d51SJames Wright 
1039ab09d51SJames Wright   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q_ceed, NULL);
1049ab09d51SJames Wright   CeedElemRestrictionCreateVector(elem_restr_grad_velo, &rhs_ceed, NULL);
1059ab09d51SJames Wright 
1069ab09d51SJames Wright   PetscCall(
1079ab09d51SJames Wright       OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, q_ceed, rhs_ceed, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
1089ab09d51SJames Wright 
1099ab09d51SJames Wright   // -- Build Mass operator
1109ab09d51SJames Wright   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
1119ab09d51SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
1129ab09d51SJames Wright   CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
1139ab09d51SJames Wright   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
1149ab09d51SJames Wright   CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
1159ab09d51SJames Wright 
1169ab09d51SJames Wright   {  // -- Setup KSP for L^2 projection with lumped mass operator
1179ab09d51SJames Wright     PetscInt l_size, g_size;
1189ab09d51SJames Wright     Mat      mat_mass;
1199ab09d51SJames Wright     VecType  vec_type;
1209ab09d51SJames Wright     Vec      M_inv;
1219ab09d51SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
1229ab09d51SJames Wright 
1239ab09d51SJames Wright     PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &M_inv));
1249ab09d51SJames Wright     PetscCall(VecGetLocalSize(M_inv, &l_size));
1259ab09d51SJames Wright     PetscCall(VecGetSize(M_inv, &g_size));
1269ab09d51SJames Wright     PetscCall(VecGetType(M_inv, &vec_type));
1279ab09d51SJames Wright     PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &M_inv));
1289ab09d51SJames Wright 
1299ab09d51SJames Wright     CeedElemRestrictionCreateVector(elem_restr_grad_velo, &mass_output, NULL);
1309ab09d51SJames Wright     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, rhs_ceed, mass_output, NULL, NULL, &mass_matop_ctx));
1319ab09d51SJames Wright     CeedVectorDestroy(&mass_output);
1329ab09d51SJames Wright 
1339ab09d51SJames Wright     PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, mass_matop_ctx, &mat_mass));
1349ab09d51SJames Wright     PetscCall(MatShellSetContextDestroy(mat_mass, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
1359ab09d51SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
1369ab09d51SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
1379ab09d51SJames Wright     PetscCall(MatShellSetVecType(mat_mass, vec_type));
1389ab09d51SJames Wright 
1399ab09d51SJames Wright     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
1409ab09d51SJames Wright     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
1419ab09d51SJames Wright     {
1429ab09d51SJames Wright       PC pc;
1439ab09d51SJames Wright       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
1449ab09d51SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
1459ab09d51SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1469ab09d51SJames Wright       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
1479ab09d51SJames Wright       // TODO Not sure if the option below are necessary
1489ab09d51SJames Wright       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
1499ab09d51SJames Wright     }
1509ab09d51SJames Wright     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
1519ab09d51SJames Wright     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
1529ab09d51SJames Wright   }
1539ab09d51SJames Wright 
1549ab09d51SJames Wright   CeedVectorDestroy(&q_ceed);
1559ab09d51SJames Wright   CeedVectorDestroy(&rhs_ceed);
1569ab09d51SJames Wright   CeedBasisDestroy(&basis_grad_velo);
1579ab09d51SJames Wright   CeedElemRestrictionDestroy(&elem_restr_grad_velo);
1589ab09d51SJames Wright   CeedQFunctionDestroy(&qf_rhs_assemble);
1599ab09d51SJames Wright   CeedQFunctionDestroy(&qf_mass);
1609ab09d51SJames Wright   CeedOperatorDestroy(&op_rhs_assemble);
1619ab09d51SJames Wright   CeedOperatorDestroy(&op_mass);
1629ab09d51SJames Wright   PetscFunctionReturn(0);
1639ab09d51SJames Wright }
164*a2401be5SJames Wright 
165*a2401be5SJames Wright PetscErrorCode VelocityGradientProjectionApply(User user, Vec Q_loc, Vec VelocityGradient) {
166*a2401be5SJames Wright   NodalProjectionData  grad_velo_proj = user->grad_velo_proj;
167*a2401be5SJames Wright   OperatorApplyContext l2_rhs_ctx     = grad_velo_proj->l2_rhs_ctx;
168*a2401be5SJames Wright 
169*a2401be5SJames Wright   PetscFunctionBeginUser;
170*a2401be5SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
171*a2401be5SJames Wright 
172*a2401be5SJames Wright   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
173*a2401be5SJames Wright 
174*a2401be5SJames Wright   PetscFunctionReturn(0);
175*a2401be5SJames Wright }
176