162b7942eSJames Wright // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors. 262b7942eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 362b7942eSJames Wright // 462b7942eSJames Wright // SPDX-License-Identifier: BSD-2-Clause 562b7942eSJames Wright // 662b7942eSJames Wright // This file is part of CEED: http://github.com/ceed 762b7942eSJames Wright 862b7942eSJames Wright #include "../qfunctions/sgs_dd_model.h" 962b7942eSJames Wright 1062b7942eSJames Wright #include <petscdmplex.h> 1162b7942eSJames Wright 1262b7942eSJames Wright #include "../navierstokes.h" 1362b7942eSJames Wright 14c38c977aSJames Wright typedef struct { 15c38c977aSJames Wright CeedElemRestriction elem_restr_grid_aniso, elem_restr_sgs; 16c38c977aSJames Wright CeedVector grid_aniso_ceed; 17*ee1455b7SJames Wright CeedQFunctionContext sgsdd_qfctx; 18c38c977aSJames Wright } *SGS_DD_ModelSetupData; 19c38c977aSJames Wright 20c38c977aSJames Wright PetscErrorCode SGS_DD_ModelSetupDataDestroy(SGS_DD_ModelSetupData sgs_dd_setup_data) { 21c38c977aSJames Wright PetscFunctionBeginUser; 22c38c977aSJames Wright CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso); 23c38c977aSJames Wright CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs); 24c38c977aSJames Wright CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed); 25*ee1455b7SJames Wright CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx); 26c38c977aSJames Wright 27c38c977aSJames Wright PetscCall(PetscFree(sgs_dd_setup_data)); 28c38c977aSJames Wright PetscFunctionReturn(0); 29c38c977aSJames Wright } 30c38c977aSJames Wright 31*ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes 32*ee1455b7SJames Wright PetscErrorCode SGS_DD_ModelCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt *num_components) { 33*ee1455b7SJames Wright PetscFE fe; 34*ee1455b7SJames Wright PetscSection section; 35*ee1455b7SJames Wright PetscInt dim; 36*ee1455b7SJames Wright 37*ee1455b7SJames Wright PetscFunctionBeginUser; 38*ee1455b7SJames Wright *num_components = 6; 39*ee1455b7SJames Wright 40*ee1455b7SJames Wright PetscCall(DMClone(dm_source, dm_sgs)); 41*ee1455b7SJames Wright PetscCall(DMGetDimension(*dm_sgs, &dim)); 42*ee1455b7SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection")); 43*ee1455b7SJames Wright 44*ee1455b7SJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, *num_components, PETSC_FALSE, degree, PETSC_DECIDE, &fe)); 45*ee1455b7SJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "Subgrid Stress Projection")); 46*ee1455b7SJames Wright PetscCall(DMAddField(*dm_sgs, NULL, (PetscObject)fe)); 47*ee1455b7SJames Wright PetscCall(DMCreateDS(*dm_sgs)); 48*ee1455b7SJames Wright PetscCall(DMPlexSetClosurePermutationTensor(*dm_sgs, PETSC_DETERMINE, NULL)); 49*ee1455b7SJames Wright 50*ee1455b7SJames Wright PetscCall(DMGetLocalSection(*dm_sgs, §ion)); 51*ee1455b7SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 52*ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX")); 53*ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY")); 54*ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ")); 55*ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ")); 56*ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ")); 57*ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY")); 58*ee1455b7SJames Wright 59*ee1455b7SJames Wright PetscCall(PetscFEDestroy(&fe)); 60*ee1455b7SJames Wright 61*ee1455b7SJames Wright PetscFunctionReturn(0); 62*ee1455b7SJames Wright }; 63*ee1455b7SJames Wright 64*ee1455b7SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes 65*ee1455b7SJames Wright PetscErrorCode SGS_DD_ModelSetupNodalEvaluation(Ceed ceed, User user, CeedData ceed_data, SGS_DD_ModelSetupData sgs_dd_setup_data) { 66*ee1455b7SJames Wright SGS_DD_Data sgs_dd_data = user->sgs_dd_data; 67*ee1455b7SJames Wright CeedQFunction qf_multiplicity, qf_sgs_dd_nodal; 68*ee1455b7SJames Wright CeedOperator op_multiplicity, op_sgs_dd_nodal; 69*ee1455b7SJames Wright CeedInt num_elem, elem_size, num_comp_q, dim, num_qpts_1d, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso; 70*ee1455b7SJames Wright CeedVector multiplicity, inv_multiplicity, grad_velo_ceed; 71*ee1455b7SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs; 72*ee1455b7SJames Wright 73*ee1455b7SJames Wright PetscFunctionBeginUser; 74*ee1455b7SJames Wright PetscCall(DMGetDimension(user->dm, &dim)); 75*ee1455b7SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x); 76*ee1455b7SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q); 77*ee1455b7SJames Wright CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso); 78*ee1455b7SJames Wright CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &num_elem); 79*ee1455b7SJames Wright CeedElemRestrictionGetElementSize(ceed_data->elem_restr_q, &elem_size); 80*ee1455b7SJames Wright CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &num_qpts_1d); 81*ee1455b7SJames Wright 82*ee1455b7SJames Wright { // Get velocity gradient information 83*ee1455b7SJames Wright CeedOperatorField op_field; 84*ee1455b7SJames Wright CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field); 85*ee1455b7SJames Wright CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo); 86*ee1455b7SJames Wright CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo); 87*ee1455b7SJames Wright CeedElemRestrictionCreateVector(elem_restr_grad_velo, &grad_velo_ceed, NULL); 88*ee1455b7SJames Wright } 89*ee1455b7SJames Wright 90*ee1455b7SJames Wright PetscCall(GetRestrictionForDomain(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, num_qpts_1d, 0, &elem_restr_sgs, NULL, NULL)); 91*ee1455b7SJames Wright CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL); 92*ee1455b7SJames Wright 93*ee1455b7SJames Wright // -- Create inverse multiplicity for correcting nodal assembly 94*ee1455b7SJames Wright CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL); 95*ee1455b7SJames Wright CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity); 96*ee1455b7SJames Wright CeedElemRestrictionCreateStrided(ceed, num_elem, elem_size, 1, num_elem * elem_size, CEED_STRIDES_BACKEND, &elem_restr_inv_multiplicity); 97*ee1455b7SJames Wright CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL); 98*ee1455b7SJames Wright 99*ee1455b7SJames Wright CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity); 100*ee1455b7SJames Wright CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE); 101*ee1455b7SJames Wright CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE); 102*ee1455b7SJames Wright 103*ee1455b7SJames Wright CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity); 104*ee1455b7SJames Wright CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling"); 105*ee1455b7SJames Wright CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 106*ee1455b7SJames Wright CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 107*ee1455b7SJames Wright CeedOperatorSetNumQuadraturePoints(op_multiplicity, elem_size); 108*ee1455b7SJames Wright 109*ee1455b7SJames Wright CeedOperatorApply(op_multiplicity, multiplicity, inv_multiplicity, CEED_REQUEST_IMMEDIATE); 110*ee1455b7SJames Wright 111*ee1455b7SJames Wright // -- Create operator for SGS DD model nodal evaluation 112*ee1455b7SJames Wright switch (user->phys->state_var) { 113*ee1455b7SJames Wright case STATEVAR_PRIMITIVE: 114*ee1455b7SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicNodal_Prim, ComputeSGS_DDAnisotropicNodal_Prim_loc, &qf_sgs_dd_nodal); 115*ee1455b7SJames Wright break; 116*ee1455b7SJames Wright case STATEVAR_CONSERVATIVE: 117*ee1455b7SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicNodal_Conserv, ComputeSGS_DDAnisotropicNodal_Conserv_loc, &qf_sgs_dd_nodal); 118*ee1455b7SJames Wright break; 119*ee1455b7SJames Wright default: 120*ee1455b7SJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, 121*ee1455b7SJames Wright "Anisotropic data-driven SGS nodal evaluation not available for chosen state variable"); 122*ee1455b7SJames Wright } 123*ee1455b7SJames Wright 124*ee1455b7SJames Wright // Mesh/geometry order and solution basis order may differ, therefore must interpolate 125*ee1455b7SJames Wright CeedBasis basis_x_to_q; 126*ee1455b7SJames Wright PetscCall(CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q)); 127*ee1455b7SJames Wright 128*ee1455b7SJames Wright CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx); 129*ee1455b7SJames Wright CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE); 130*ee1455b7SJames Wright CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP); 131*ee1455b7SJames Wright CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE); 132*ee1455b7SJames Wright CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE); 133*ee1455b7SJames Wright CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE); 134*ee1455b7SJames Wright CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE); 135*ee1455b7SJames Wright 136*ee1455b7SJames Wright CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal); 137*ee1455b7SJames Wright CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, user->q_ceed); 138*ee1455b7SJames Wright CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord); 139*ee1455b7SJames Wright CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 140*ee1455b7SJames Wright CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_COLLOCATED, 141*ee1455b7SJames Wright sgs_dd_setup_data->grid_aniso_ceed); 142*ee1455b7SJames Wright CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, inv_multiplicity); 143*ee1455b7SJames Wright CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 144*ee1455b7SJames Wright 145*ee1455b7SJames Wright PetscCall(OperatorApplyContextCreate(user->grad_velo_proj->dm, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_nodal, grad_velo_ceed, 146*ee1455b7SJames Wright sgs_dd_data->sgs_nodal_ceed, NULL, NULL, &sgs_dd_data->op_nodal_evaluation_ctx)); 147*ee1455b7SJames Wright 148*ee1455b7SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 149*ee1455b7SJames Wright 150*ee1455b7SJames Wright CeedVectorDestroy(&multiplicity); 151*ee1455b7SJames Wright CeedVectorDestroy(&inv_multiplicity); 152*ee1455b7SJames Wright CeedVectorDestroy(&grad_velo_ceed); 153*ee1455b7SJames Wright CeedBasisDestroy(&basis_x_to_q); 154*ee1455b7SJames Wright CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity); 155*ee1455b7SJames Wright CeedQFunctionDestroy(&qf_multiplicity); 156*ee1455b7SJames Wright CeedQFunctionDestroy(&qf_sgs_dd_nodal); 157*ee1455b7SJames Wright CeedOperatorDestroy(&op_multiplicity); 158*ee1455b7SJames Wright CeedOperatorDestroy(&op_sgs_dd_nodal); 159*ee1455b7SJames Wright PetscFunctionReturn(0); 160*ee1455b7SJames Wright } 161*ee1455b7SJames Wright 16262b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN 16362b7942eSJames Wright PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) { 16462b7942eSJames Wright PetscFunctionBeginUser; 16562b7942eSJames Wright for (PetscInt i = 0; i < N; i++) { 16662b7942eSJames Wright for (PetscInt j = 0; j < M; j++) { 16762b7942eSJames Wright B[j * N + i] = A[i * M + j]; 16862b7942eSJames Wright } 16962b7942eSJames Wright } 17062b7942eSJames Wright PetscFunctionReturn(0); 17162b7942eSJames Wright } 17262b7942eSJames Wright 17362b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct 17462b7942eSJames Wright PetscErrorCode SGS_DD_ModelContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SGS_DDModelContext *psgsdd_ctx) { 17562b7942eSJames Wright SGS_DDModelContext sgsdd_ctx; 17662b7942eSJames Wright PetscInt num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons; 17762b7942eSJames Wright char file_path[PETSC_MAX_PATH_LEN]; 17862b7942eSJames Wright PetscScalar *temp; 17962b7942eSJames Wright 18062b7942eSJames Wright PetscFunctionBeginUser; 18162b7942eSJames Wright { 18262b7942eSJames Wright SGS_DDModelContext sgsdd_temp; 18362b7942eSJames Wright PetscCall(PetscNew(&sgsdd_temp)); 18462b7942eSJames Wright *sgsdd_temp = **psgsdd_ctx; 18562b7942eSJames Wright sgsdd_temp->offsets.bias1 = 0; 18662b7942eSJames Wright sgsdd_temp->offsets.bias2 = sgsdd_temp->offsets.bias1 + num_neurons; 18762b7942eSJames Wright sgsdd_temp->offsets.weight1 = sgsdd_temp->offsets.bias2 + num_neurons; 18862b7942eSJames Wright sgsdd_temp->offsets.weight2 = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs; 18962b7942eSJames Wright sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons; 19062b7942eSJames Wright PetscInt total_num_scalars = sgsdd_temp->offsets.out_scaling + 2 * num_outputs; 19162b7942eSJames Wright sgsdd_temp->total_bytes = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]); 19262b7942eSJames Wright PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx)); 19362b7942eSJames Wright *sgsdd_ctx = *sgsdd_temp; 19462b7942eSJames Wright PetscCall(PetscFree(sgsdd_temp)); 19562b7942eSJames Wright } 19662b7942eSJames Wright 19762b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat")); 19862b7942eSJames Wright PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1])); 19962b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat")); 20062b7942eSJames Wright PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2])); 20162b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat")); 20262b7942eSJames Wright PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling])); 20362b7942eSJames Wright 20462b7942eSJames Wright { 20562b7942eSJames Wright PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp)); 20662b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat")); 20762b7942eSJames Wright PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp)); 20862b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons)); 20962b7942eSJames Wright PetscCall(PetscFree(temp)); 21062b7942eSJames Wright } 21162b7942eSJames Wright { 21262b7942eSJames Wright PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp)); 21362b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat")); 21462b7942eSJames Wright PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp)); 21562b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs)); 21662b7942eSJames Wright PetscCall(PetscFree(temp)); 21762b7942eSJames Wright } 21862b7942eSJames Wright 21962b7942eSJames Wright PetscCall(PetscFree(*psgsdd_ctx)); 22062b7942eSJames Wright *psgsdd_ctx = sgsdd_ctx; 22162b7942eSJames Wright PetscFunctionReturn(0); 22262b7942eSJames Wright } 22362b7942eSJames Wright 22462b7942eSJames Wright PetscErrorCode SGS_DD_ModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 225*ee1455b7SJames Wright PetscReal alpha = 0; 22662b7942eSJames Wright SGS_DDModelContext sgsdd_ctx; 22762b7942eSJames Wright MPI_Comm comm = user->comm; 228*ee1455b7SJames Wright char sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters"; 229c38c977aSJames Wright SGS_DD_ModelSetupData sgs_dd_setup_data; 230*ee1455b7SJames Wright NewtonianIdealGasContext gas; 23162b7942eSJames Wright PetscFunctionBeginUser; 23262b7942eSJames Wright 2339ab09d51SJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem)); 2349ab09d51SJames Wright 23562b7942eSJames Wright PetscCall(PetscNew(&sgsdd_ctx)); 23662b7942eSJames Wright 23762b7942eSJames Wright PetscOptionsBegin(comm, NULL, "SGS Data-Drive Model Options", NULL); 23862b7942eSJames Wright PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL)); 23962b7942eSJames Wright PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir, 24062b7942eSJames Wright sgs_dd_dir, sizeof(sgs_dd_dir), NULL)); 24162b7942eSJames Wright PetscOptionsEnd(); 24262b7942eSJames Wright 24362b7942eSJames Wright sgsdd_ctx->num_layers = 2; 24462b7942eSJames Wright sgsdd_ctx->num_inputs = 6; 24562b7942eSJames Wright sgsdd_ctx->num_outputs = 6; 24662b7942eSJames Wright sgsdd_ctx->num_neurons = 20; 24762b7942eSJames Wright sgsdd_ctx->alpha = alpha; 24862b7942eSJames Wright 24962b7942eSJames Wright // PetscCall(SGS_DD_ModelContextFill(comm, sgs_dd_dir, &sgsdd_ctx)); 25062b7942eSJames Wright 251*ee1455b7SJames Wright // -- Create DM for storing SGS tensor at nodes 252*ee1455b7SJames Wright PetscCall(PetscNew(&user->sgs_dd_data)); 253*ee1455b7SJames Wright PetscCall(SGS_DD_ModelCreateDM(user->dm, &user->sgs_dd_data->dm_sgs, user->app_ctx->degree, &user->sgs_dd_data->num_comp_sgs)); 254*ee1455b7SJames Wright 255c38c977aSJames Wright PetscCall(PetscNew(&sgs_dd_setup_data)); 256c38c977aSJames Wright 257*ee1455b7SJames Wright CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas); 258*ee1455b7SJames Wright sgsdd_ctx->gas = *gas; 259*ee1455b7SJames Wright CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas); 260*ee1455b7SJames Wright CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx); 261*ee1455b7SJames Wright CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx); 262*ee1455b7SJames Wright CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc); 263*ee1455b7SJames Wright 264c38c977aSJames Wright // -- Compute and store anisotropy tensor 265c38c977aSJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso, 266c38c977aSJames Wright &sgs_dd_setup_data->grid_aniso_ceed)); 267c38c977aSJames Wright 268*ee1455b7SJames Wright // -- Create Nodal Evaluation Operator 269*ee1455b7SJames Wright PetscCall(SGS_DD_ModelSetupNodalEvaluation(ceed, user, ceed_data, sgs_dd_setup_data)); 270*ee1455b7SJames Wright 271c38c977aSJames Wright PetscCall(SGS_DD_ModelSetupDataDestroy(sgs_dd_setup_data)); 27262b7942eSJames Wright PetscFunctionReturn(0); 27362b7942eSJames Wright } 274*ee1455b7SJames Wright 275*ee1455b7SJames Wright PetscErrorCode SGS_DD_DataDestroy(SGS_DD_Data sgs_dd_data) { 276*ee1455b7SJames Wright PetscFunctionBeginUser; 277*ee1455b7SJames Wright if (!sgs_dd_data) PetscFunctionReturn(0); 278*ee1455b7SJames Wright 279*ee1455b7SJames Wright CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed); 280*ee1455b7SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx)); 281*ee1455b7SJames Wright PetscCall(DMDestroy(&sgs_dd_data->dm_sgs)); 282*ee1455b7SJames Wright PetscCall(PetscFree(sgs_dd_data)); 283*ee1455b7SJames Wright 284*ee1455b7SJames Wright PetscFunctionReturn(0); 285*ee1455b7SJames Wright } 286