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