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