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