1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2530ad8c4SKenneth E. Jansen // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3530ad8c4SKenneth E. Jansen // 4530ad8c4SKenneth E. Jansen // SPDX-License-Identifier: BSD-2-Clause 5530ad8c4SKenneth E. Jansen // 6530ad8c4SKenneth E. Jansen // This file is part of CEED: http://github.com/ceed 7530ad8c4SKenneth E. Jansen 8530ad8c4SKenneth E. Jansen /// @file 9530ad8c4SKenneth E. Jansen /// Utility functions for setting up Gaussian Wave problem 10530ad8c4SKenneth E. Jansen 11530ad8c4SKenneth E. Jansen #include "../qfunctions/gaussianwave.h" 12530ad8c4SKenneth E. Jansen 1349aac155SJeremy L Thompson #include <ceed.h> 1449aac155SJeremy L Thompson #include <petscdm.h> 1549aac155SJeremy L Thompson 16530ad8c4SKenneth E. Jansen #include "../navierstokes.h" 179f844368SJames Wright #include "../qfunctions/bc_freestream_type.h" 18530ad8c4SKenneth E. Jansen 19530ad8c4SKenneth E. Jansen PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 20530ad8c4SKenneth E. Jansen User user = *(User *)ctx; 21a424bcd0SJames Wright MPI_Comm comm = user->comm; 22a424bcd0SJames Wright Ceed ceed = user->ceed; 23530ad8c4SKenneth E. Jansen GaussianWaveContext gausswave_ctx; 24530ad8c4SKenneth E. Jansen FreestreamContext freestream_ctx; 25530ad8c4SKenneth E. Jansen NewtonianIdealGasContext newtonian_ig_ctx; 26530ad8c4SKenneth E. Jansen CeedQFunctionContext gausswave_context; 27530ad8c4SKenneth E. Jansen 28530ad8c4SKenneth E. Jansen PetscFunctionBeginUser; 29530ad8c4SKenneth E. Jansen PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 30530ad8c4SKenneth E. Jansen 31530ad8c4SKenneth E. Jansen switch (user->phys->state_var) { 32530ad8c4SKenneth E. Jansen case STATEVAR_CONSERVATIVE: 33530ad8c4SKenneth E. Jansen problem->ics.qfunction = IC_GaussianWave_Conserv; 34530ad8c4SKenneth E. Jansen problem->ics.qfunction_loc = IC_GaussianWave_Conserv_loc; 35530ad8c4SKenneth E. Jansen break; 36530ad8c4SKenneth E. Jansen case STATEVAR_PRIMITIVE: 37530ad8c4SKenneth E. Jansen problem->ics.qfunction = IC_GaussianWave_Prim; 38530ad8c4SKenneth E. Jansen problem->ics.qfunction_loc = IC_GaussianWave_Prim_loc; 39530ad8c4SKenneth E. Jansen break; 40530ad8c4SKenneth E. Jansen } 41530ad8c4SKenneth E. Jansen 42530ad8c4SKenneth E. Jansen // -- Option Defaults 43530ad8c4SKenneth E. Jansen CeedScalar epicenter[3] = {0.}; // m 44530ad8c4SKenneth E. Jansen CeedScalar width = 0.002; // m 45530ad8c4SKenneth E. Jansen CeedScalar amplitude = 0.1; // - 46530ad8c4SKenneth E. Jansen 47530ad8c4SKenneth E. Jansen PetscOptionsBegin(comm, NULL, "Options for GAUSSIAN_WAVE problem", NULL); 48530ad8c4SKenneth E. Jansen PetscInt narray = 3; 49530ad8c4SKenneth E. Jansen PetscCall(PetscOptionsScalarArray("-epicenter", "Coordinates of center of wave", NULL, epicenter, &narray, NULL)); 50530ad8c4SKenneth E. Jansen PetscCheck(narray == 3, comm, PETSC_ERR_ARG_SIZ, "-epicenter should recieve array of size 3, instead recieved size %" PetscInt_FMT ".", narray); 51530ad8c4SKenneth E. Jansen PetscCall(PetscOptionsScalar("-width", "Width parameter for perturbation size", NULL, width, &width, NULL)); 52530ad8c4SKenneth E. Jansen PetscCall(PetscOptionsScalar("-amplitude", "Amplitude of the perturbation", NULL, amplitude, &litude, NULL)); 53530ad8c4SKenneth E. Jansen PetscOptionsEnd(); 54530ad8c4SKenneth E. Jansen 55530ad8c4SKenneth E. Jansen width *= user->units->meter; 56530ad8c4SKenneth E. Jansen for (int i = 0; i < 3; i++) epicenter[i] *= user->units->meter; 57530ad8c4SKenneth E. Jansen 58530ad8c4SKenneth E. Jansen // -- Set gausswave_ctx struct values 59530ad8c4SKenneth E. Jansen PetscCall(PetscCalloc1(1, &gausswave_ctx)); 60a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 61a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_freestream.qfunction_context, CEED_MEM_HOST, &freestream_ctx)); 62530ad8c4SKenneth E. Jansen 63530ad8c4SKenneth E. Jansen gausswave_ctx->amplitude = amplitude; 64530ad8c4SKenneth E. Jansen gausswave_ctx->width = width; 65530ad8c4SKenneth E. Jansen gausswave_ctx->S_infty = freestream_ctx->S_infty; 66530ad8c4SKenneth E. Jansen gausswave_ctx->newt_ctx = *newtonian_ig_ctx; 67530ad8c4SKenneth E. Jansen PetscCall(PetscArraycpy(gausswave_ctx->epicenter, epicenter, 3)); 68530ad8c4SKenneth E. Jansen 69a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 70a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_freestream.qfunction_context, &freestream_ctx)); 71530ad8c4SKenneth E. Jansen 72a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &gausswave_context)); 73a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(gausswave_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx)); 74a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(gausswave_context, CEED_MEM_HOST, FreeContextPetsc)); 75a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 76530ad8c4SKenneth E. Jansen problem->ics.qfunction_context = gausswave_context; 77ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 78530ad8c4SKenneth E. Jansen } 79