1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Utility functions for setting up Gaussian Wave problem 6 7 #include "../qfunctions/gaussianwave.h" 8 9 #include <ceed.h> 10 #include <petscdm.h> 11 12 #include <navierstokes.h> 13 #include "../qfunctions/bc_freestream_type.h" 14 15 PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx) { 16 Honee honee = *(Honee *)ctx; 17 MPI_Comm comm = honee->comm; 18 Ceed ceed = honee->ceed; 19 GaussianWaveContext gausswave_ctx; 20 NewtonianIdealGasContext newtonian_ig_ctx; 21 CeedQFunctionContext gausswave_qfctx; 22 23 PetscFunctionBeginUser; 24 PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx)); 25 26 switch (honee->phys->state_var) { 27 case STATEVAR_CONSERVATIVE: 28 problem->ics.qf_func_ptr = IC_GaussianWave_Conserv; 29 problem->ics.qf_loc = IC_GaussianWave_Conserv_loc; 30 break; 31 case STATEVAR_PRIMITIVE: 32 problem->ics.qf_func_ptr = IC_GaussianWave_Prim; 33 problem->ics.qf_loc = IC_GaussianWave_Prim_loc; 34 break; 35 case STATEVAR_ENTROPY: 36 problem->ics.qf_func_ptr = IC_GaussianWave_Entropy; 37 problem->ics.qf_loc = IC_GaussianWave_Entropy_loc; 38 break; 39 } 40 41 // -- Option Defaults 42 CeedScalar epicenter[3] = {0.}; // m 43 CeedScalar width = 0.002; // m 44 CeedScalar amplitude = 0.1; // - 45 46 PetscOptionsBegin(comm, NULL, "Options for GAUSSIAN_WAVE problem", NULL); 47 PetscInt narray = 3; 48 PetscCall(PetscOptionsScalarArray("-epicenter", "Coordinates of center of wave", NULL, epicenter, &narray, NULL)); 49 PetscCheck(narray == 3, comm, PETSC_ERR_ARG_SIZ, "-epicenter should recieve array of size 3, instead recieved size %" PetscInt_FMT ".", narray); 50 PetscCall(PetscOptionsScalar("-width", "Width parameter for perturbation size", NULL, width, &width, NULL)); 51 PetscCall(PetscOptionsScalar("-amplitude", "Amplitude of the perturbation", NULL, amplitude, &litude, NULL)); 52 PetscOptionsEnd(); 53 54 width *= honee->units->meter; 55 for (int i = 0; i < 3; i++) epicenter[i] *= honee->units->meter; 56 57 // -- Set gausswave_ctx struct values 58 PetscCall(PetscNew(&gausswave_ctx)); 59 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 60 for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 61 BCDefinition bc_def = problem->bc_defs[b]; 62 HoneeBCStruct honee_bc; 63 const char *name; 64 65 PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 66 if (!strcmp(name, "freestream")) { 67 FreestreamContext freestream_ctx; 68 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 69 PetscCallCeed(ceed, CeedQFunctionContextGetData(honee_bc->qfctx, CEED_MEM_HOST, &freestream_ctx)); 70 gausswave_ctx->S_infty = freestream_ctx->S_infty; 71 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(honee_bc->qfctx, &freestream_ctx)); 72 break; 73 } 74 } 75 76 gausswave_ctx->amplitude = amplitude; 77 gausswave_ctx->width = width; 78 gausswave_ctx->newt_ctx = *newtonian_ig_ctx; 79 PetscCall(PetscArraycpy(gausswave_ctx->epicenter, epicenter, 3)); 80 81 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 82 83 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &gausswave_qfctx)); 84 PetscCallCeed(ceed, CeedQFunctionContextSetData(gausswave_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx)); 85 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(gausswave_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 86 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 87 problem->ics.qfctx = gausswave_qfctx; 88 PetscFunctionReturn(PETSC_SUCCESS); 89 } 90