1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 #include "../qfunctions/grid_anisotropy_tensor.h" 5 6 #include <petscdmplex.h> 7 8 #include <navierstokes.h> 9 10 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 11 CeedVector *grid_aniso_vector) { 12 NodalProjectionData grid_aniso_proj; 13 CeedBasis basis_grid_aniso; 14 CeedVector q_data; 15 CeedElemRestriction elem_restr_qd; 16 CeedInt q_data_size; 17 MPI_Comm comm = PetscObjectComm((PetscObject)honee->dm); 18 PetscInt height = 0, dm_field = 0; 19 20 PetscFunctionBeginUser; 21 PetscCall(PetscNew(&grid_aniso_proj)); 22 23 { // -- Create DM for Anisotropic tensor L^2 projection 24 PetscSection section; 25 26 PetscCall(DMClone(honee->dm, &grid_aniso_proj->dm)); 27 PetscCall(DMSetMatrixPreallocateSkip(grid_aniso_proj->dm, PETSC_TRUE)); 28 PetscCall(PetscObjectSetName((PetscObject)grid_aniso_proj->dm, "Grid Anisotropy Tensor Projection")); 29 30 // -- Setup DM 31 grid_aniso_proj->num_comp = 7; 32 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, 1, &grid_aniso_proj->num_comp, 33 grid_aniso_proj->dm)); 34 35 PetscCall(DMGetLocalSection(grid_aniso_proj->dm, §ion)); 36 PetscCall(PetscSectionSetFieldName(section, 0, "")); 37 PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMGridAnisotropyTensorXX")); 38 PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMGridAnisotropyTensorYY")); 39 PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMGridAnisotropyTensorZZ")); 40 PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMGridAnisotropyTensorYZ")); 41 PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMGridAnisotropyTensorXZ")); 42 PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMGridAnisotropyTensorXY")); 43 PetscCall(PetscSectionSetComponentName(section, 0, 6, "GridAnisotropyTensorFrobNorm")); 44 } 45 46 // -- Get Pre-requisite things 47 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grid_aniso_proj->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, 48 elem_restr_grid_aniso)); 49 PetscCall(DMPlexCeedBasisCreate(ceed, grid_aniso_proj->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_grid_aniso)); 50 PetscCall(QDataGet(ceed, grid_aniso_proj->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 51 &elem_restr_qd, &q_data, &q_data_size)); 52 53 { // -- Build RHS operator 54 CeedOperator op_rhs_assemble; 55 CeedQFunction qf_rhs_assemble; 56 57 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, AnisotropyTensorProjection, AnisotropyTensorProjection_loc, &qf_rhs_assemble)); 58 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE)); 59 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "v", grid_aniso_proj->num_comp, CEED_EVAL_INTERP)); 60 61 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble)); 62 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 63 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "v", *elem_restr_grid_aniso, basis_grid_aniso, CEED_VECTOR_ACTIVE)); 64 65 PetscCall(OperatorApplyContextCreate(honee->dm, grid_aniso_proj->dm, ceed, op_rhs_assemble, CEED_VECTOR_NONE, NULL, NULL, NULL, 66 &grid_aniso_proj->l2_rhs_ctx)); 67 68 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble)); 69 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble)); 70 } 71 72 { // Setup KSP for L^2 projection 73 CeedOperator op_mass; 74 CeedQFunction qf_mass; 75 Mat mat_mass; 76 77 PetscCall(HoneeMassQFunctionCreate(ceed, grid_aniso_proj->num_comp, q_data_size, &qf_mass)); 78 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 79 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", *elem_restr_grid_aniso, basis_grid_aniso, CEED_VECTOR_ACTIVE)); 80 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 81 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", *elem_restr_grid_aniso, basis_grid_aniso, CEED_VECTOR_ACTIVE)); 82 83 PetscCall(MatCreateCeed(grid_aniso_proj->dm, grid_aniso_proj->dm, op_mass, NULL, &mat_mass)); 84 85 PetscCall(KSPCreate(comm, &grid_aniso_proj->ksp)); 86 PetscCall(KSPSetOptionsPrefix(grid_aniso_proj->ksp, "grid_anisotropy_tensor_projection_")); 87 { 88 PC pc; 89 PetscCall(KSPGetPC(grid_aniso_proj->ksp, &pc)); 90 PetscCall(PCSetType(pc, PCJACOBI)); 91 PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 92 PetscCall(KSPSetType(grid_aniso_proj->ksp, KSPCG)); 93 PetscCall(KSPSetNormType(grid_aniso_proj->ksp, KSP_NORM_NATURAL)); 94 PetscCall(KSPSetTolerances(grid_aniso_proj->ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 95 } 96 PetscCall(KSPSetFromOptions_WithMatCeed(grid_aniso_proj->ksp, mat_mass)); 97 98 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 99 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 100 PetscCall(MatDestroy(&mat_mass)); 101 } 102 103 { // -- Project anisotropy data and store in CeedVector 104 Vec Grid_Anisotropy, grid_anisotropy_loc; 105 106 PetscCall(DMGetGlobalVector(grid_aniso_proj->dm, &Grid_Anisotropy)); 107 PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Grid_Anisotropy, grid_aniso_proj->l2_rhs_ctx)); 108 PetscCall(KSPSolve(grid_aniso_proj->ksp, Grid_Anisotropy, Grid_Anisotropy)); 109 110 // Copy anisotropy tensor data to CeedVector 111 PetscCall(DMGetLocalVector(grid_aniso_proj->dm, &grid_anisotropy_loc)); 112 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(*elem_restr_grid_aniso, grid_aniso_vector, NULL)); 113 PetscCall(DMGlobalToLocal(grid_aniso_proj->dm, Grid_Anisotropy, INSERT_VALUES, grid_anisotropy_loc)); 114 PetscCall(VecCopyPetscToCeed(grid_anisotropy_loc, *grid_aniso_vector)); 115 PetscCall(DMRestoreLocalVector(grid_aniso_proj->dm, &grid_anisotropy_loc)); 116 PetscCall(DMRestoreGlobalVector(grid_aniso_proj->dm, &Grid_Anisotropy)); 117 } 118 119 PetscCall(NodalProjectionDataDestroy(&grid_aniso_proj)); 120 PetscCallCeed(ceed, CeedBasisDestroy(&basis_grid_aniso)); 121 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 122 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 123 PetscFunctionReturn(PETSC_SUCCESS); 124 } 125 126 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso, 127 CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso) { 128 CeedInt q_data_size, num_nodes; 129 CeedQFunction qf_colloc; 130 CeedOperator op_colloc; 131 CeedVector q_data; 132 CeedElemRestriction elem_restr_qd; 133 PetscInt height = 0; 134 135 PetscFunctionBeginUser; 136 *num_comp_aniso = 7; 137 { 138 CeedBasis basis_q; 139 PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q)); 140 PetscCallCeed(ceed, CeedBasisGetNumNodes(basis_q, &num_nodes)); 141 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 142 } 143 PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, *num_comp_aniso, 144 elem_restr_grid_aniso)); 145 PetscCall(QDataGet(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, 146 &q_data, &q_data_size)); 147 148 // -- Build collocation operator 149 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, AnisotropyTensorCollocate, AnisotropyTensorCollocate_loc, &qf_colloc)); 150 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_colloc, "qdata", q_data_size, CEED_EVAL_NONE)); 151 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_colloc, "v", *num_comp_aniso, CEED_EVAL_NONE)); 152 153 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_colloc, NULL, NULL, &op_colloc)); 154 PetscCallCeed(ceed, CeedOperatorSetField(op_colloc, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 155 PetscCallCeed(ceed, CeedOperatorSetField(op_colloc, "v", *elem_restr_grid_aniso, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 156 157 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(*elem_restr_grid_aniso, aniso_colloc_ceed, NULL)); 158 159 PetscCallCeed(ceed, CeedOperatorApply(op_colloc, CEED_VECTOR_NONE, *aniso_colloc_ceed, CEED_REQUEST_IMMEDIATE)); 160 161 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 162 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 163 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_colloc)); 164 PetscCallCeed(ceed, CeedOperatorDestroy(&op_colloc)); 165 PetscFunctionReturn(PETSC_SUCCESS); 166 } 167