1*530ad8c4SKenneth E. Jansen // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2*530ad8c4SKenneth E. Jansen // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3*530ad8c4SKenneth E. Jansen // 4*530ad8c4SKenneth E. Jansen // SPDX-License-Identifier: BSD-2-Clause 5*530ad8c4SKenneth E. Jansen // 6*530ad8c4SKenneth E. Jansen // This file is part of CEED: http://github.com/ceed 7*530ad8c4SKenneth E. Jansen 8*530ad8c4SKenneth E. Jansen /// @file 9*530ad8c4SKenneth E. Jansen /// Utility functions for setting up Gaussian Wave problem 10*530ad8c4SKenneth E. Jansen 11*530ad8c4SKenneth E. Jansen #include "../qfunctions/gaussianwave.h" 12*530ad8c4SKenneth E. Jansen 13*530ad8c4SKenneth E. Jansen #include "../navierstokes.h" 14*530ad8c4SKenneth E. Jansen #include "../qfunctions/freestream_bc_type.h" 15*530ad8c4SKenneth E. Jansen 16*530ad8c4SKenneth E. Jansen PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 17*530ad8c4SKenneth E. Jansen User user = *(User *)ctx; 18*530ad8c4SKenneth E. Jansen MPI_Comm comm = PETSC_COMM_WORLD; 19*530ad8c4SKenneth E. Jansen GaussianWaveContext gausswave_ctx; 20*530ad8c4SKenneth E. Jansen FreestreamContext freestream_ctx; 21*530ad8c4SKenneth E. Jansen NewtonianIdealGasContext newtonian_ig_ctx; 22*530ad8c4SKenneth E. Jansen CeedQFunctionContext gausswave_context; 23*530ad8c4SKenneth E. Jansen 24*530ad8c4SKenneth E. Jansen PetscFunctionBeginUser; 25*530ad8c4SKenneth E. Jansen PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 26*530ad8c4SKenneth E. Jansen 27*530ad8c4SKenneth E. Jansen switch (user->phys->state_var) { 28*530ad8c4SKenneth E. Jansen case STATEVAR_CONSERVATIVE: 29*530ad8c4SKenneth E. Jansen problem->ics.qfunction = IC_GaussianWave_Conserv; 30*530ad8c4SKenneth E. Jansen problem->ics.qfunction_loc = IC_GaussianWave_Conserv_loc; 31*530ad8c4SKenneth E. Jansen break; 32*530ad8c4SKenneth E. Jansen case STATEVAR_PRIMITIVE: 33*530ad8c4SKenneth E. Jansen problem->ics.qfunction = IC_GaussianWave_Prim; 34*530ad8c4SKenneth E. Jansen problem->ics.qfunction_loc = IC_GaussianWave_Prim_loc; 35*530ad8c4SKenneth E. Jansen break; 36*530ad8c4SKenneth E. Jansen } 37*530ad8c4SKenneth E. Jansen 38*530ad8c4SKenneth E. Jansen // -- Option Defaults 39*530ad8c4SKenneth E. Jansen CeedScalar epicenter[3] = {0.}; // m 40*530ad8c4SKenneth E. Jansen CeedScalar width = 0.002; // m 41*530ad8c4SKenneth E. Jansen CeedScalar amplitude = 0.1; // - 42*530ad8c4SKenneth E. Jansen 43*530ad8c4SKenneth E. Jansen PetscOptionsBegin(comm, NULL, "Options for GAUSSIAN_WAVE problem", NULL); 44*530ad8c4SKenneth E. Jansen PetscInt narray = 3; 45*530ad8c4SKenneth E. Jansen PetscCall(PetscOptionsScalarArray("-epicenter", "Coordinates of center of wave", NULL, epicenter, &narray, NULL)); 46*530ad8c4SKenneth E. Jansen PetscCheck(narray == 3, comm, PETSC_ERR_ARG_SIZ, "-epicenter should recieve array of size 3, instead recieved size %" PetscInt_FMT ".", narray); 47*530ad8c4SKenneth E. Jansen PetscCall(PetscOptionsScalar("-width", "Width parameter for perturbation size", NULL, width, &width, NULL)); 48*530ad8c4SKenneth E. Jansen PetscCall(PetscOptionsScalar("-amplitude", "Amplitude of the perturbation", NULL, amplitude, &litude, NULL)); 49*530ad8c4SKenneth E. Jansen PetscOptionsEnd(); 50*530ad8c4SKenneth E. Jansen 51*530ad8c4SKenneth E. Jansen width *= user->units->meter; 52*530ad8c4SKenneth E. Jansen for (int i = 0; i < 3; i++) epicenter[i] *= user->units->meter; 53*530ad8c4SKenneth E. Jansen 54*530ad8c4SKenneth E. Jansen // -- Set gausswave_ctx struct values 55*530ad8c4SKenneth E. Jansen PetscCall(PetscCalloc1(1, &gausswave_ctx)); 56*530ad8c4SKenneth E. Jansen CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 57*530ad8c4SKenneth E. Jansen CeedQFunctionContextGetData(problem->apply_freestream.qfunction_context, CEED_MEM_HOST, &freestream_ctx); 58*530ad8c4SKenneth E. Jansen 59*530ad8c4SKenneth E. Jansen gausswave_ctx->amplitude = amplitude; 60*530ad8c4SKenneth E. Jansen gausswave_ctx->width = width; 61*530ad8c4SKenneth E. Jansen gausswave_ctx->S_infty = freestream_ctx->S_infty; 62*530ad8c4SKenneth E. Jansen gausswave_ctx->newt_ctx = *newtonian_ig_ctx; 63*530ad8c4SKenneth E. Jansen PetscCall(PetscArraycpy(gausswave_ctx->epicenter, epicenter, 3)); 64*530ad8c4SKenneth E. Jansen 65*530ad8c4SKenneth E. Jansen CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 66*530ad8c4SKenneth E. Jansen CeedQFunctionContextRestoreData(problem->apply_freestream.qfunction_context, &freestream_ctx); 67*530ad8c4SKenneth E. Jansen 68*530ad8c4SKenneth E. Jansen CeedQFunctionContextCreate(user->ceed, &gausswave_context); 69*530ad8c4SKenneth E. Jansen CeedQFunctionContextSetData(gausswave_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx); 70*530ad8c4SKenneth E. Jansen CeedQFunctionContextSetDataDestroy(gausswave_context, CEED_MEM_HOST, FreeContextPetsc); 71*530ad8c4SKenneth E. Jansen CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 72*530ad8c4SKenneth E. Jansen problem->ics.qfunction_context = gausswave_context; 73*530ad8c4SKenneth E. Jansen 74*530ad8c4SKenneth E. Jansen PetscFunctionReturn(0); 75*530ad8c4SKenneth E. Jansen } 76