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 26e7754af5SKenneth E. Jansen // -- Option Defaults 27e7754af5SKenneth E. Jansen CeedScalar epicenter[3] = {0.}; // m 28e7754af5SKenneth E. Jansen CeedScalar width = 0.002; // m 29e7754af5SKenneth E. Jansen CeedScalar amplitude = 0.1; // - 30e7754af5SKenneth E. Jansen 31e7754af5SKenneth E. Jansen PetscOptionsBegin(comm, NULL, "Options for GAUSSIAN_WAVE problem", NULL); 32e7754af5SKenneth E. Jansen PetscInt narray = 3; 33e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalarArray("-epicenter", "Coordinates of center of wave", NULL, epicenter, &narray, NULL)); 34e7754af5SKenneth E. Jansen PetscCheck(narray == 3, comm, PETSC_ERR_ARG_SIZ, "-epicenter should recieve array of size 3, instead recieved size %" PetscInt_FMT ".", narray); 35e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-width", "Width parameter for perturbation size", NULL, width, &width, NULL)); 36e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-amplitude", "Amplitude of the perturbation", NULL, amplitude, &litude, NULL)); 37e7754af5SKenneth E. Jansen PetscOptionsEnd(); 38e7754af5SKenneth E. Jansen 390c373b74SJames Wright width *= honee->units->meter; 400c373b74SJames Wright for (int i = 0; i < 3; i++) epicenter[i] *= honee->units->meter; 41e7754af5SKenneth E. Jansen 42e7754af5SKenneth E. Jansen // -- Set gausswave_ctx struct values 432d898fa6SJames Wright PetscCall(PetscNew(&gausswave_ctx)); 44e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 45d4e423e7SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 46d4e423e7SJames Wright BCDefinition bc_def = problem->bc_defs[b]; 47d4e423e7SJames Wright HoneeBCStruct honee_bc; 48d4e423e7SJames Wright const char *name; 49d4e423e7SJames Wright 50d4e423e7SJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 51d4e423e7SJames Wright if (!strcmp(name, "freestream")) { 52d4e423e7SJames Wright FreestreamContext freestream_ctx; 53d4e423e7SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 54d4e423e7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(honee_bc->qfctx, CEED_MEM_HOST, &freestream_ctx)); 55d4e423e7SJames Wright gausswave_ctx->S_infty = freestream_ctx->S_infty; 56d4e423e7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(honee_bc->qfctx, &freestream_ctx)); 57d4e423e7SJames Wright break; 58d4e423e7SJames Wright } 59d4e423e7SJames Wright } 60e7754af5SKenneth E. Jansen 61e7754af5SKenneth E. Jansen gausswave_ctx->amplitude = amplitude; 62e7754af5SKenneth E. Jansen gausswave_ctx->width = width; 63e7754af5SKenneth E. Jansen gausswave_ctx->newt_ctx = *newtonian_ig_ctx; 64e7754af5SKenneth E. Jansen PetscCall(PetscArraycpy(gausswave_ctx->epicenter, epicenter, 3)); 65e7754af5SKenneth E. Jansen 66e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 67e7754af5SKenneth E. Jansen 680c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &gausswave_qfctx)); 69e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(gausswave_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx)); 70e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(gausswave_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 71e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 72*f5dc303cSJames Wright 73*f5dc303cSJames Wright switch (honee->phys->state_var) { 74*f5dc303cSJames Wright case STATEVAR_CONSERVATIVE: 75*f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = IC_GaussianWave_Conserv, .qf_loc = IC_GaussianWave_Conserv_loc}; 76*f5dc303cSJames Wright break; 77*f5dc303cSJames Wright case STATEVAR_PRIMITIVE: 78*f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = IC_GaussianWave_Prim, .qf_loc = IC_GaussianWave_Prim_loc}; 79*f5dc303cSJames Wright break; 80*f5dc303cSJames Wright case STATEVAR_ENTROPY: 81*f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = IC_GaussianWave_Entropy, .qf_loc = IC_GaussianWave_Entropy_loc}; 82*f5dc303cSJames Wright break; 83*f5dc303cSJames Wright } 84e07531f7SJames Wright problem->ics.qfctx = gausswave_qfctx; 85d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 86e7754af5SKenneth E. Jansen } 87