xref: /libCEED/examples/fluids/src/velocity_gradient_projection.c (revision 8e6aa226c2c84e58dd7feb551fd506c4f25986db)
1 // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 /// @file
8 /// Functions for setting up and projecting the velocity gradient
9 
10 #include "../qfunctions/velocity_gradient_projection.h"
11 
12 #include <petscdmplex.h>
13 
14 #include "../navierstokes.h"
15 
16 PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) {
17   PetscFE      fe;
18   PetscSection section;
19   PetscInt     dim;
20 
21   PetscFunctionBeginUser;
22   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
23 
24   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
25   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
26   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
27 
28   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, grad_velo_proj->num_comp, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
29   PetscCall(PetscObjectSetName((PetscObject)fe, "Velocity Gradient Projection"));
30   PetscCall(DMAddField(grad_velo_proj->dm, NULL, (PetscObject)fe));
31   PetscCall(DMCreateDS(grad_velo_proj->dm));
32   PetscCall(DMPlexSetClosurePermutationTensor(grad_velo_proj->dm, PETSC_DETERMINE, NULL));
33 
34   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
35   PetscCall(PetscSectionSetFieldName(section, 0, ""));
36   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
37   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
38   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
39   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
40   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
41   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
42   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
43   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
44   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
45 
46   PetscCall(PetscFEDestroy(&fe));
47   PetscFunctionReturn(PETSC_SUCCESS);
48 };
49 
50 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
51   OperatorApplyContext mass_matop_ctx;
52   CeedOperator         op_rhs_assemble, op_mass;
53   CeedQFunction        qf_rhs_assemble, qf_mass;
54   CeedBasis            basis_grad_velo;
55   CeedElemRestriction  elem_restr_grad_velo;
56   PetscInt             dim, num_comp_x, num_comp_q, q_data_size, num_qpts_1d, num_nodes_1d;
57 
58   PetscFunctionBeginUser;
59   PetscCall(PetscNew(&user->grad_velo_proj));
60   NodalProjectionData grad_velo_proj = user->grad_velo_proj;
61 
62   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
63 
64   // -- Get Pre-requisite things
65   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
66   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
67   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
68   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &num_qpts_1d);
69   CeedBasisGetNumNodes1D(ceed_data->basis_q, &num_nodes_1d);
70   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
71 
72   PetscCall(GetRestrictionForDomain(ceed, grad_velo_proj->dm, 0, 0, 0, 0, num_qpts_1d, q_data_size, &elem_restr_grad_velo, NULL, NULL));
73 
74   CeedBasisCreateTensorH1Lagrange(ceed, dim, grad_velo_proj->num_comp, num_nodes_1d, num_qpts_1d, CEED_GAUSS, &basis_grad_velo);
75 
76   // -- Build RHS operator
77   switch (user->phys->state_var) {
78     case STATEVAR_PRIMITIVE:
79       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble);
80       break;
81     case STATEVAR_CONSERVATIVE:
82       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, &qf_rhs_assemble);
83       break;
84     default:
85       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
86   }
87 
88   CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context);
89   CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_q, CEED_EVAL_INTERP);
90   CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
91   CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE);
92   CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP);
93   CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP);
94 
95   CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble);
96   CeedOperatorSetField(op_rhs_assemble, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
97   CeedOperatorSetField(op_rhs_assemble, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
98   CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
99   CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
100   CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
101 
102   PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
103 
104   // -- Build Mass operator
105   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
106   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
107   CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
108   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
109   CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
110 
111   {  // -- Setup KSP for L^2 projection with lumped mass operator
112     Mat      mat_mass;
113     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
114 
115     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, NULL, NULL, NULL, NULL, &mass_matop_ctx));
116     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
117 
118     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
119     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
120     {
121       PC pc;
122       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
123       PetscCall(PCSetType(pc, PCJACOBI));
124       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
125       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
126       // TODO Not sure if the option below are necessary
127       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
128     }
129     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
130     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
131   }
132 
133   CeedBasisDestroy(&basis_grad_velo);
134   CeedElemRestrictionDestroy(&elem_restr_grad_velo);
135   CeedQFunctionDestroy(&qf_rhs_assemble);
136   CeedQFunctionDestroy(&qf_mass);
137   CeedOperatorDestroy(&op_rhs_assemble);
138   CeedOperatorDestroy(&op_mass);
139   PetscFunctionReturn(PETSC_SUCCESS);
140 }
141 
142 PetscErrorCode VelocityGradientProjectionApply(User user, Vec Q_loc, Vec VelocityGradient) {
143   NodalProjectionData  grad_velo_proj = user->grad_velo_proj;
144   OperatorApplyContext l2_rhs_ctx     = grad_velo_proj->l2_rhs_ctx;
145 
146   PetscFunctionBeginUser;
147   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
148 
149   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
150 
151   PetscFunctionReturn(PETSC_SUCCESS);
152 }
153