xref: /honee/problems/gaussianwave.c (revision d949ddfc92e3b3e8aad90700e3fb60ca7f4585db)
1e7754af5SKenneth E. Jansen // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2e7754af5SKenneth E. Jansen // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3e7754af5SKenneth E. Jansen //
4e7754af5SKenneth E. Jansen // SPDX-License-Identifier: BSD-2-Clause
5e7754af5SKenneth E. Jansen //
6e7754af5SKenneth E. Jansen // This file is part of CEED:  http://github.com/ceed
7e7754af5SKenneth E. Jansen 
8e7754af5SKenneth E. Jansen /// @file
9e7754af5SKenneth E. Jansen /// Utility functions for setting up Gaussian Wave problem
10e7754af5SKenneth E. Jansen 
11e7754af5SKenneth E. Jansen #include "../qfunctions/gaussianwave.h"
12e7754af5SKenneth E. Jansen 
13e419654dSJeremy L Thompson #include <ceed.h>
14e419654dSJeremy L Thompson #include <petscdm.h>
15e419654dSJeremy L Thompson 
16e7754af5SKenneth E. Jansen #include "../navierstokes.h"
17e7754af5SKenneth E. Jansen #include "../qfunctions/freestream_bc_type.h"
18e7754af5SKenneth E. Jansen 
19e7754af5SKenneth E. Jansen PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
20e7754af5SKenneth E. Jansen   User                     user = *(User *)ctx;
21e7754af5SKenneth E. Jansen   MPI_Comm                 comm = PETSC_COMM_WORLD;
22e7754af5SKenneth E. Jansen   GaussianWaveContext      gausswave_ctx;
23e7754af5SKenneth E. Jansen   FreestreamContext        freestream_ctx;
24e7754af5SKenneth E. Jansen   NewtonianIdealGasContext newtonian_ig_ctx;
25e7754af5SKenneth E. Jansen   CeedQFunctionContext     gausswave_context;
26e7754af5SKenneth E. Jansen 
27e7754af5SKenneth E. Jansen   PetscFunctionBeginUser;
28e7754af5SKenneth E. Jansen   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
29e7754af5SKenneth E. Jansen 
30e7754af5SKenneth E. Jansen   switch (user->phys->state_var) {
31e7754af5SKenneth E. Jansen     case STATEVAR_CONSERVATIVE:
32e7754af5SKenneth E. Jansen       problem->ics.qfunction     = IC_GaussianWave_Conserv;
33e7754af5SKenneth E. Jansen       problem->ics.qfunction_loc = IC_GaussianWave_Conserv_loc;
34e7754af5SKenneth E. Jansen       break;
35e7754af5SKenneth E. Jansen     case STATEVAR_PRIMITIVE:
36e7754af5SKenneth E. Jansen       problem->ics.qfunction     = IC_GaussianWave_Prim;
37e7754af5SKenneth E. Jansen       problem->ics.qfunction_loc = IC_GaussianWave_Prim_loc;
38e7754af5SKenneth E. Jansen       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, &amplitude, NULL));
52e7754af5SKenneth E. Jansen   PetscOptionsEnd();
53e7754af5SKenneth E. Jansen 
54e7754af5SKenneth E. Jansen   width *= user->units->meter;
55e7754af5SKenneth E. Jansen   for (int i = 0; i < 3; i++) epicenter[i] *= user->units->meter;
56e7754af5SKenneth E. Jansen 
57e7754af5SKenneth E. Jansen   // -- Set gausswave_ctx struct values
58e7754af5SKenneth E. Jansen   PetscCall(PetscCalloc1(1, &gausswave_ctx));
59e7754af5SKenneth E. Jansen   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
60e7754af5SKenneth E. Jansen   CeedQFunctionContextGetData(problem->apply_freestream.qfunction_context, CEED_MEM_HOST, &freestream_ctx);
61e7754af5SKenneth E. Jansen 
62e7754af5SKenneth E. Jansen   gausswave_ctx->amplitude = amplitude;
63e7754af5SKenneth E. Jansen   gausswave_ctx->width     = width;
64e7754af5SKenneth E. Jansen   gausswave_ctx->S_infty   = freestream_ctx->S_infty;
65e7754af5SKenneth E. Jansen   gausswave_ctx->newt_ctx  = *newtonian_ig_ctx;
66e7754af5SKenneth E. Jansen   PetscCall(PetscArraycpy(gausswave_ctx->epicenter, epicenter, 3));
67e7754af5SKenneth E. Jansen 
68e7754af5SKenneth E. Jansen   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
69e7754af5SKenneth E. Jansen   CeedQFunctionContextRestoreData(problem->apply_freestream.qfunction_context, &freestream_ctx);
70e7754af5SKenneth E. Jansen 
71e7754af5SKenneth E. Jansen   CeedQFunctionContextCreate(user->ceed, &gausswave_context);
72e7754af5SKenneth E. Jansen   CeedQFunctionContextSetData(gausswave_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx);
73e7754af5SKenneth E. Jansen   CeedQFunctionContextSetDataDestroy(gausswave_context, CEED_MEM_HOST, FreeContextPetsc);
74e7754af5SKenneth E. Jansen   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
75e7754af5SKenneth E. Jansen   problem->ics.qfunction_context = gausswave_context;
76e7754af5SKenneth E. Jansen 
77*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
78e7754af5SKenneth E. Jansen }
79