xref: /libCEED/examples/fluids/src/velocity_gradient_projection.c (revision f84dcd568de24bead0416446b432c2186c5e82c6)
1999ff5c7SJames Wright // Copyright (c) 2017-2023, 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   PetscFE      fe;
18999ff5c7SJames Wright   PetscSection section;
19999ff5c7SJames Wright   PetscInt     dim;
20999ff5c7SJames Wright 
21999ff5c7SJames Wright   PetscFunctionBeginUser;
22999ff5c7SJames Wright   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
23999ff5c7SJames Wright 
24999ff5c7SJames Wright   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
25999ff5c7SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
26999ff5c7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
27999ff5c7SJames Wright 
28999ff5c7SJames Wright   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, grad_velo_proj->num_comp, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
29999ff5c7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "Velocity Gradient Projection"));
30999ff5c7SJames Wright   PetscCall(DMAddField(grad_velo_proj->dm, NULL, (PetscObject)fe));
31999ff5c7SJames Wright   PetscCall(DMCreateDS(grad_velo_proj->dm));
32999ff5c7SJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(grad_velo_proj->dm, PETSC_DETERMINE, NULL));
33999ff5c7SJames Wright 
34999ff5c7SJames Wright   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
35999ff5c7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
36999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
37999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
38999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
39999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
40999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
41999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
42999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
43999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
44999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
45999ff5c7SJames Wright 
46999ff5c7SJames Wright   PetscCall(PetscFEDestroy(&fe));
47999ff5c7SJames Wright   PetscFunctionReturn(0);
48999ff5c7SJames Wright };
493909d146SJames Wright 
503909d146SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
513909d146SJames Wright   OperatorApplyContext mass_matop_ctx;
523909d146SJames Wright   CeedOperator         op_rhs_assemble, op_mass;
533909d146SJames Wright   CeedQFunction        qf_rhs_assemble, qf_mass;
543909d146SJames Wright   CeedBasis            basis_grad_velo;
553909d146SJames Wright   CeedVector           q_ceed, rhs_ceed, mass_output;
563909d146SJames Wright   CeedElemRestriction  elem_restr_grad_velo;
573909d146SJames Wright   PetscInt             dim, num_comp_x, num_comp_q, q_data_size, num_qpts_1d, num_nodes_1d;
583909d146SJames Wright 
593909d146SJames Wright   PetscFunctionBeginUser;
603909d146SJames Wright   PetscCall(PetscNew(&user->grad_velo_proj));
613909d146SJames Wright   NodalProjectionData grad_velo_proj = user->grad_velo_proj;
623909d146SJames Wright 
633909d146SJames Wright   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
643909d146SJames Wright 
653909d146SJames Wright   // -- Get Pre-requisite things
663909d146SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
673909d146SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
683909d146SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
693909d146SJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &num_qpts_1d);
703909d146SJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &num_nodes_1d);
713909d146SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
723909d146SJames Wright 
733909d146SJames Wright   PetscCall(GetRestrictionForDomain(ceed, grad_velo_proj->dm, 0, 0, 0, num_qpts_1d, q_data_size, &elem_restr_grad_velo, NULL, NULL));
743909d146SJames Wright 
753909d146SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, grad_velo_proj->num_comp, num_nodes_1d, num_qpts_1d, CEED_GAUSS, &basis_grad_velo);
763909d146SJames Wright 
773909d146SJames Wright   // -- Build RHS operator
783909d146SJames Wright   switch (user->phys->state_var) {
793909d146SJames Wright     case STATEVAR_PRIMITIVE:
803909d146SJames Wright       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble);
813909d146SJames Wright       break;
823909d146SJames Wright     case STATEVAR_CONSERVATIVE:
833909d146SJames Wright       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, &qf_rhs_assemble);
843909d146SJames Wright       break;
853909d146SJames Wright     default:
863909d146SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
873909d146SJames Wright   }
883909d146SJames Wright 
893909d146SJames Wright   CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context);
903909d146SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_q, CEED_EVAL_INTERP);
913909d146SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
923909d146SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE);
933909d146SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP);
943909d146SJames Wright   CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP);
953909d146SJames Wright 
963909d146SJames Wright   CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble);
973909d146SJames Wright   CeedOperatorSetField(op_rhs_assemble, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
983909d146SJames Wright   CeedOperatorSetField(op_rhs_assemble, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
993909d146SJames Wright   CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
1003909d146SJames Wright   CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
1013909d146SJames Wright   CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
1023909d146SJames Wright 
1033909d146SJames Wright   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q_ceed, NULL);
1043909d146SJames Wright   CeedElemRestrictionCreateVector(elem_restr_grad_velo, &rhs_ceed, NULL);
1053909d146SJames Wright 
1063909d146SJames Wright   PetscCall(
1073909d146SJames Wright       OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, q_ceed, rhs_ceed, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
1083909d146SJames Wright 
1093909d146SJames Wright   // -- Build Mass operator
1103909d146SJames Wright   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
1113909d146SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
1123909d146SJames Wright   CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
1133909d146SJames Wright   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
1143909d146SJames Wright   CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
1153909d146SJames Wright 
1163909d146SJames Wright   {  // -- Setup KSP for L^2 projection with lumped mass operator
1173909d146SJames Wright     PetscInt l_size, g_size;
1183909d146SJames Wright     Mat      mat_mass;
1193909d146SJames Wright     VecType  vec_type;
1203909d146SJames Wright     Vec      M_inv;
1213909d146SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
1223909d146SJames Wright 
1233909d146SJames Wright     PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &M_inv));
1243909d146SJames Wright     PetscCall(VecGetLocalSize(M_inv, &l_size));
1253909d146SJames Wright     PetscCall(VecGetSize(M_inv, &g_size));
1263909d146SJames Wright     PetscCall(VecGetType(M_inv, &vec_type));
1273909d146SJames Wright     PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &M_inv));
1283909d146SJames Wright 
1293909d146SJames Wright     CeedElemRestrictionCreateVector(elem_restr_grad_velo, &mass_output, NULL);
1303909d146SJames Wright     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, rhs_ceed, mass_output, NULL, NULL, &mass_matop_ctx));
1313909d146SJames Wright     CeedVectorDestroy(&mass_output);
1323909d146SJames Wright 
1333909d146SJames Wright     PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, mass_matop_ctx, &mat_mass));
1343909d146SJames Wright     PetscCall(MatShellSetContextDestroy(mat_mass, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
1353909d146SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
1363909d146SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
1373909d146SJames Wright     PetscCall(MatShellSetVecType(mat_mass, vec_type));
1383909d146SJames Wright 
1393909d146SJames Wright     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
1403909d146SJames Wright     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
1413909d146SJames Wright     {
1423909d146SJames Wright       PC pc;
1433909d146SJames Wright       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
1443909d146SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
1453909d146SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1463909d146SJames Wright       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
1473909d146SJames Wright       // TODO Not sure if the option below are necessary
1483909d146SJames Wright       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
1493909d146SJames Wright     }
1503909d146SJames Wright     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
1513909d146SJames Wright     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
1523909d146SJames Wright   }
1533909d146SJames Wright 
1543909d146SJames Wright   CeedVectorDestroy(&q_ceed);
1553909d146SJames Wright   CeedVectorDestroy(&rhs_ceed);
1563909d146SJames Wright   CeedBasisDestroy(&basis_grad_velo);
1573909d146SJames Wright   CeedElemRestrictionDestroy(&elem_restr_grad_velo);
1583909d146SJames Wright   CeedQFunctionDestroy(&qf_rhs_assemble);
1593909d146SJames Wright   CeedQFunctionDestroy(&qf_mass);
1603909d146SJames Wright   CeedOperatorDestroy(&op_rhs_assemble);
1613909d146SJames Wright   CeedOperatorDestroy(&op_mass);
1623909d146SJames Wright   PetscFunctionReturn(0);
1633909d146SJames Wright }
164*f84dcd56SJames Wright 
165*f84dcd56SJames Wright PetscErrorCode VelocityGradientProjectionApply(User user, Vec Q_loc, Vec VelocityGradient) {
166*f84dcd56SJames Wright   NodalProjectionData  grad_velo_proj = user->grad_velo_proj;
167*f84dcd56SJames Wright   OperatorApplyContext l2_rhs_ctx     = grad_velo_proj->l2_rhs_ctx;
168*f84dcd56SJames Wright 
169*f84dcd56SJames Wright   PetscFunctionBeginUser;
170*f84dcd56SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
171*f84dcd56SJames Wright 
172*f84dcd56SJames Wright   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
173*f84dcd56SJames Wright 
174*f84dcd56SJames Wright   PetscFunctionReturn(0);
175*f84dcd56SJames Wright }
176