xref: /libCEED/examples/fluids/navierstokes.c (revision 7659d40c2b62cf54bbe42ee118e449c0e8801b81)
1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4ccaff030SJeremy L Thompson //
5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9ccaff030SJeremy L Thompson //
10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16ccaff030SJeremy L Thompson 
17ccaff030SJeremy L Thompson //                        libCEED + PETSc Example: Navier-Stokes
18ccaff030SJeremy L Thompson //
19ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve a
20ccaff030SJeremy L Thompson // Navier-Stokes problem.
21ccaff030SJeremy L Thompson //
22ccaff030SJeremy L Thompson // The code is intentionally "raw", using only low-level communication
23ccaff030SJeremy L Thompson // primitives.
24ccaff030SJeremy L Thompson //
25ccaff030SJeremy L Thompson // Build with:
26ccaff030SJeremy L Thompson //
27ccaff030SJeremy L Thompson //     make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes
28ccaff030SJeremy L Thompson //
29ccaff030SJeremy L Thompson // Sample runs:
30ccaff030SJeremy L Thompson //
3184d34d69SLeila Ghaffari //     ./navierstokes -ceed /cpu/self -problem density_current -degree 1
3284d34d69SLeila Ghaffari //     ./navierstokes -ceed /gpu/occa -problem advection -degree 1
33ccaff030SJeremy L Thompson //
3484d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test explicit -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3
3584d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test implicit_stab_none -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha
3684d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test implicit_stab_supg -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -stab supg
37ccaff030SJeremy L Thompson 
38ccaff030SJeremy L Thompson /// @file
39ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc
40ccaff030SJeremy L Thompson 
41ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n";
42ccaff030SJeremy L Thompson 
43ccaff030SJeremy L Thompson #include <petscts.h>
44ccaff030SJeremy L Thompson #include <petscdmplex.h>
45ccaff030SJeremy L Thompson #include <ceed.h>
46ccaff030SJeremy L Thompson #include <stdbool.h>
47ccaff030SJeremy L Thompson #include <petscsys.h>
48ccaff030SJeremy L Thompson #include "common.h"
49b0137797SLeila Ghaffari #include "setup-boundary.h"
50ccaff030SJeremy L Thompson #include "advection.h"
51ccaff030SJeremy L Thompson #include "advection2d.h"
52ccaff030SJeremy L Thompson #include "densitycurrent.h"
53ccaff030SJeremy L Thompson 
5484d34d69SLeila Ghaffari #if PETSC_VERSION_LT(3,14,0)
5584d34d69SLeila Ghaffari #  define DMPlexGetClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexGetClosureIndices(a,b,c,d,f,g,i)
5684d34d69SLeila Ghaffari #  define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i)
5784d34d69SLeila Ghaffari #endif
5884d34d69SLeila Ghaffari 
5984d34d69SLeila Ghaffari // MemType Options
6084d34d69SLeila Ghaffari static const char *const memTypes[] = {
6184d34d69SLeila Ghaffari   "host",
6284d34d69SLeila Ghaffari   "device",
6384d34d69SLeila Ghaffari   "memType", "CEED_MEM_", NULL
6484d34d69SLeila Ghaffari };
6584d34d69SLeila Ghaffari 
66ccaff030SJeremy L Thompson // Problem Options
67ccaff030SJeremy L Thompson typedef enum {
68ccaff030SJeremy L Thompson   NS_DENSITY_CURRENT = 0,
69ccaff030SJeremy L Thompson   NS_ADVECTION = 1,
70ccaff030SJeremy L Thompson   NS_ADVECTION2D = 2,
71ccaff030SJeremy L Thompson } problemType;
72ccaff030SJeremy L Thompson static const char *const problemTypes[] = {
73ccaff030SJeremy L Thompson   "density_current",
74ccaff030SJeremy L Thompson   "advection",
75ccaff030SJeremy L Thompson   "advection2d",
7684d34d69SLeila Ghaffari   "problemType", "NS_", NULL
77ccaff030SJeremy L Thompson };
78ccaff030SJeremy L Thompson 
791184866aSLeila Ghaffari // Wind Options for Advection
801184866aSLeila Ghaffari typedef enum {
811184866aSLeila Ghaffari   ADVECTION_WIND_ROTATION = 0,
821184866aSLeila Ghaffari   ADVECTION_WIND_TRANSLATION = 1,
831184866aSLeila Ghaffari } WindType;
841184866aSLeila Ghaffari static const char *const WindTypes[] = {
851184866aSLeila Ghaffari   "rotation",
861184866aSLeila Ghaffari   "translation",
871184866aSLeila Ghaffari   "WindType", "ADVECTION_WIND_", NULL
881184866aSLeila Ghaffari };
891184866aSLeila Ghaffari 
90ccaff030SJeremy L Thompson typedef enum {
91ccaff030SJeremy L Thompson   STAB_NONE = 0,
92ccaff030SJeremy L Thompson   STAB_SU = 1,   // Streamline Upwind
93ccaff030SJeremy L Thompson   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
94ccaff030SJeremy L Thompson } StabilizationType;
95ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = {
9684d34d69SLeila Ghaffari   "none",
97ccaff030SJeremy L Thompson   "SU",
98ccaff030SJeremy L Thompson   "SUPG",
99ccaff030SJeremy L Thompson   "StabilizationType", "STAB_", NULL
100ccaff030SJeremy L Thompson };
101ccaff030SJeremy L Thompson 
10284d34d69SLeila Ghaffari // Test Options
10384d34d69SLeila Ghaffari typedef enum {
10484d34d69SLeila Ghaffari   TEST_NONE = 0,               // Non test mode
10584d34d69SLeila Ghaffari   TEST_EXPLICIT = 1,           // Explicit test
10684d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_NONE = 2, // Implicit test no stab
10784d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_SUPG = 3, // Implicit test supg stab
10884d34d69SLeila Ghaffari } testType;
10984d34d69SLeila Ghaffari static const char *const testTypes[] = {
11084d34d69SLeila Ghaffari   "none",
11184d34d69SLeila Ghaffari   "explicit",
11284d34d69SLeila Ghaffari   "implicit_stab_none",
11384d34d69SLeila Ghaffari   "implicit_stab_supg",
11484d34d69SLeila Ghaffari   "testType", "TEST_", NULL
11584d34d69SLeila Ghaffari };
11684d34d69SLeila Ghaffari 
11784d34d69SLeila Ghaffari // Tests specific data
11884d34d69SLeila Ghaffari typedef struct {
11984d34d69SLeila Ghaffari   PetscScalar testtol;
12084d34d69SLeila Ghaffari   const char *filepath;
12184d34d69SLeila Ghaffari } testData;
12284d34d69SLeila Ghaffari 
12384d34d69SLeila Ghaffari testData testOptions[] = {
12484d34d69SLeila Ghaffari   [TEST_NONE] = {
12584d34d69SLeila Ghaffari     .testtol = 0.,
12684d34d69SLeila Ghaffari     .filepath = NULL
12784d34d69SLeila Ghaffari   },
12884d34d69SLeila Ghaffari   [TEST_EXPLICIT] = {
12984d34d69SLeila Ghaffari     .testtol = 1E-5,
13084d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-explicit.bin"
13184d34d69SLeila Ghaffari   },
13284d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_NONE] = {
13384d34d69SLeila Ghaffari     .testtol = 5E-4,
13484d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-none.bin"
13584d34d69SLeila Ghaffari   },
13684d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_SUPG] = {
13784d34d69SLeila Ghaffari     .testtol = 5E-4,
13884d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-supg.bin"
13984d34d69SLeila Ghaffari   }
14084d34d69SLeila Ghaffari };
14184d34d69SLeila Ghaffari 
142ccaff030SJeremy L Thompson // Problem specific data
143ccaff030SJeremy L Thompson typedef struct {
1448b982baeSLeila Ghaffari   CeedInt dim, qdatasizeVol, qdatasizeSur;
145*7659d40cSLeila Ghaffari   CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applyVol_ifunction, applySur;
146ccaff030SJeremy L Thompson   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
147ccaff030SJeremy L Thompson                        PetscScalar[], void *);
148*7659d40cSLeila Ghaffari   const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc,
149*7659d40cSLeila Ghaffari              *applyVol_ifunction_loc, *applySur_loc;
150ccaff030SJeremy L Thompson   const bool non_zero_time;
151ccaff030SJeremy L Thompson } problemData;
152ccaff030SJeremy L Thompson 
153ccaff030SJeremy L Thompson problemData problemOptions[] = {
154ccaff030SJeremy L Thompson   [NS_DENSITY_CURRENT] = {
155ccaff030SJeremy L Thompson     .dim                       = 3,
156ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1578b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
158b0137797SLeila Ghaffari     .setupVol                  = Setup,
159b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
160356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
161356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
162ccaff030SJeremy L Thompson     .ics                       = ICsDC,
163ccaff030SJeremy L Thompson     .ics_loc                   = ICsDC_loc,
164c96c872fSLeila Ghaffari     .applyVol_rhs              = DC,
165c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = DC_loc,
166c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_DC,
167c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_DC_loc,
168ccaff030SJeremy L Thompson     .bc                        = Exact_DC,
16984d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
170ccaff030SJeremy L Thompson   },
171ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
172ccaff030SJeremy L Thompson     .dim                       = 3,
173ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1748b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
175b0137797SLeila Ghaffari     .setupVol                  = Setup,
176b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
177356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
178356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
179ccaff030SJeremy L Thompson     .ics                       = ICsAdvection,
180ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection_loc,
181c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection,
182c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection_loc,
183c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection,
184c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection_loc,
185*7659d40cSLeila Ghaffari     .applySur                  = Advection_Sur,
186*7659d40cSLeila Ghaffari     .applySur_loc              = Advection_Sur_loc,
187ccaff030SJeremy L Thompson     .bc                        = Exact_Advection,
18884d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
189ccaff030SJeremy L Thompson   },
190ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
191ccaff030SJeremy L Thompson     .dim                       = 2,
192ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 5,
1938b982baeSLeila Ghaffari     .qdatasizeSur              = 3,
194c96c872fSLeila Ghaffari     .setupVol                  = Setup2d,
195c96c872fSLeila Ghaffari     .setupVol_loc              = Setup2d_loc,
196b0137797SLeila Ghaffari     .setupSur                  = SetupBoundary2d,
197b0137797SLeila Ghaffari     .setupSur_loc              = SetupBoundary2d_loc,
198ccaff030SJeremy L Thompson     .ics                       = ICsAdvection2d,
199ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection2d_loc,
200c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection2d,
201c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection2d_loc,
202c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection2d,
203c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection2d_loc,
204*7659d40cSLeila Ghaffari     .applySur                  = Advection2d_Sur,
205*7659d40cSLeila Ghaffari     .applySur_loc              = Advection2d_Sur_loc,
206ccaff030SJeremy L Thompson     .bc                        = Exact_Advection2d,
20784d34d69SLeila Ghaffari     .non_zero_time             = PETSC_TRUE,
208ccaff030SJeremy L Thompson   },
209ccaff030SJeremy L Thompson };
210ccaff030SJeremy L Thompson 
211ccaff030SJeremy L Thompson // PETSc user data
212ccaff030SJeremy L Thompson typedef struct User_ *User;
213ccaff030SJeremy L Thompson typedef struct Units_ *Units;
214ccaff030SJeremy L Thompson 
215ccaff030SJeremy L Thompson struct User_ {
216ccaff030SJeremy L Thompson   MPI_Comm comm;
217ccaff030SJeremy L Thompson   PetscInt outputfreq;
218ccaff030SJeremy L Thompson   DM dm;
219ccaff030SJeremy L Thompson   DM dmviz;
220ccaff030SJeremy L Thompson   Mat interpviz;
221ccaff030SJeremy L Thompson   Ceed ceed;
222ccaff030SJeremy L Thompson   Units units;
223ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
2241e150236SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction;
225ccaff030SJeremy L Thompson   Vec M;
226ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
227ccaff030SJeremy L Thompson   PetscInt contsteps;
228ccaff030SJeremy L Thompson };
229ccaff030SJeremy L Thompson 
230ccaff030SJeremy L Thompson struct Units_ {
231ccaff030SJeremy L Thompson   // fundamental units
232ccaff030SJeremy L Thompson   PetscScalar meter;
233ccaff030SJeremy L Thompson   PetscScalar kilogram;
234ccaff030SJeremy L Thompson   PetscScalar second;
235ccaff030SJeremy L Thompson   PetscScalar Kelvin;
236ccaff030SJeremy L Thompson   // derived units
237ccaff030SJeremy L Thompson   PetscScalar Pascal;
238ccaff030SJeremy L Thompson   PetscScalar JperkgK;
239ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
240ccaff030SJeremy L Thompson   PetscScalar WpermK;
241ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
242ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
243ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
244ccaff030SJeremy L Thompson };
245ccaff030SJeremy L Thompson 
246ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
247ccaff030SJeremy L Thompson struct SimpleBC_ {
248*7659d40cSLeila Ghaffari   PetscInt nwall, nslip[3];
249*7659d40cSLeila Ghaffari   PetscInt walls[6], slips[3][6];
25084d34d69SLeila Ghaffari   PetscBool userbc;
251ccaff030SJeremy L Thompson };
252ccaff030SJeremy L Thompson 
253ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
254ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
255ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
256ccaff030SJeremy L Thompson }
257ccaff030SJeremy L Thompson 
258ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
259ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
26084d34d69SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value,
261ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
262ccaff030SJeremy L Thompson 
263ccaff030SJeremy L Thompson   PetscSection section;
2641184866aSLeila Ghaffari   PetscInt p, Nelem, Ndof, *erestrict, eoffset, nfields, dim, depth;
2650c6c0b13SLeila Ghaffari   DMLabel depthLabel;
2660c6c0b13SLeila Ghaffari   IS depthIS, iterIS;
26784d34d69SLeila Ghaffari   Vec Uloc;
2680c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
269ccaff030SJeremy L Thompson   PetscErrorCode ierr;
270ccaff030SJeremy L Thompson 
271ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
272ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
273da51bdd9SLeila Ghaffari   dim -= height;
274ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
275ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
276ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
277ccaff030SJeremy L Thompson   fieldoff[0] = 0;
278ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
279ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
280ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
281ccaff030SJeremy L Thompson   }
282ccaff030SJeremy L Thompson 
2830c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2840c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2850c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2860c6c0b13SLeila Ghaffari   if (domainLabel) {
2870c6c0b13SLeila Ghaffari     IS domainIS;
2880c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2890c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2900c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2910c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2920c6c0b13SLeila Ghaffari   } else {
2930c6c0b13SLeila Ghaffari     iterIS = depthIS;
2940c6c0b13SLeila Ghaffari   }
2950c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2960c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
297ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
2980c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
2990c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
300ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
30184d34d69SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
30284d34d69SLeila Ghaffari                                    &numindices, &indices, NULL, NULL);
30384d34d69SLeila Ghaffari     CHKERRQ(ierr);
30432b5ec5fSJed Brown     bool flip = false;
30532b5ec5fSJed Brown     if (height > 0) {
30632b5ec5fSJed Brown       PetscInt numCells, numFaces, start = -1;
30732b5ec5fSJed Brown       const PetscInt *orients, *faces, *cells;
30832b5ec5fSJed Brown       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
30932b5ec5fSJed Brown       ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr);
31032b5ec5fSJed Brown       if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Expected one cell in support of exterior face, but got %D cells", numCells);
31132b5ec5fSJed Brown       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
31232b5ec5fSJed Brown       ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr);
31332b5ec5fSJed Brown       for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;}
31432b5ec5fSJed Brown       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %D in cone of its support", c);
31532b5ec5fSJed Brown       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
31632b5ec5fSJed Brown       if (orients[start] < 0) flip = true;
31732b5ec5fSJed Brown     }
31884d34d69SLeila Ghaffari     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
31984d34d69SLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
32084d34d69SLeila Ghaffari           c);
321ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
322ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
32332b5ec5fSJed Brown       PetscInt ii = i;
32432b5ec5fSJed Brown       if (flip) {
32532b5ec5fSJed Brown     	  if (P == nnodes) ii = nnodes - 1 - i;
32632b5ec5fSJed Brown         else if (P*P == nnodes) {
32732b5ec5fSJed Brown           PetscInt row = i / P, col = i % P;
32832b5ec5fSJed Brown 	        ii = row + col * P;
32932b5ec5fSJed Brown 	      } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for flipping point with %D nodes != P (%D) or P^2", nnodes, P);
33032b5ec5fSJed Brown       }
331ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
332ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
333ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
334ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
33532b5ec5fSJed Brown           if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j])
33632b5ec5fSJed Brown               != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j)
337ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
338ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
33932b5ec5fSJed Brown                      c, ii, f, j);
340ccaff030SJeremy L Thompson         }
341ccaff030SJeremy L Thompson       }
342ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
34332b5ec5fSJed Brown       PetscInt loc = Involute(indices[ii*ncomp[0]]);
3446f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
345ccaff030SJeremy L Thompson     }
34684d34d69SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
34784d34d69SLeila Ghaffari                                        &numindices, &indices, NULL, NULL);
34884d34d69SLeila Ghaffari     CHKERRQ(ierr);
349ccaff030SJeremy L Thompson   }
3500c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
3510c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3520c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
353ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
3540c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
3550c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
3560c6c0b13SLeila Ghaffari 
357ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
358ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
359ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
3606f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
3616f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
3626f55dfd5Svaleriabarra                             Erestrict);
363ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
364ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
365ccaff030SJeremy L Thompson }
366ccaff030SJeremy L Thompson 
367c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
3681e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
3691e150236SLeila Ghaffari     DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
3701e150236SLeila Ghaffari     CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
3711e150236SLeila Ghaffari     CeedElemRestriction *restrictqdi) {
372c96c872fSLeila Ghaffari 
373c96c872fSLeila Ghaffari   DM dmcoord;
3741e150236SLeila Ghaffari   CeedInt dim, localNelem;
3751e150236SLeila Ghaffari   CeedInt Qdim;
376c96c872fSLeila Ghaffari   PetscErrorCode ierr;
377c96c872fSLeila Ghaffari 
378c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
3791e150236SLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
3801e150236SLeila Ghaffari   dim -= height;
3811e150236SLeila Ghaffari   Qdim = CeedIntPow(Q, dim);
382c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
383c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
384c96c872fSLeila Ghaffari   CHKERRQ(ierr);
385c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
386c96c872fSLeila Ghaffari   CHKERRQ(ierr);
387c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
388c96c872fSLeila Ghaffari   CHKERRQ(ierr);
389c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
390c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
391c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
392c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
393c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
394c96c872fSLeila Ghaffari }
395c96c872fSLeila Ghaffari 
3961e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
397*7659d40cSLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc,
398*7659d40cSLeila Ghaffari     WindType wind_type, CeedOperator op_applyVol, CeedQFunction qf_applySur,
399*7659d40cSLeila Ghaffari     CeedQFunction qf_setupSur, CeedInt height, CeedInt numP_Sur, CeedInt numQ_Sur,
400*7659d40cSLeila Ghaffari     CeedInt qdatasizeSur, CeedInt NqptsSur, CeedBasis basisxSur,
401*7659d40cSLeila Ghaffari     CeedBasis basisqSur, CeedOperator *op_apply) {
402ca3ac6ddSLeila Ghaffari 
403*7659d40cSLeila Ghaffari   CeedInt dim, nFace;
404*7659d40cSLeila Ghaffari   PetscInt lsize, localNelemSur[6];
4051e150236SLeila Ghaffari   Vec Xloc;
406*7659d40cSLeila Ghaffari   CeedVector xcorners, qdataSur[6];
407*7659d40cSLeila Ghaffari   CeedOperator op_setupSur[6], op_applySur[6];
408*7659d40cSLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
409ca3ac6ddSLeila Ghaffari   DMLabel domainLabel;
4101e150236SLeila Ghaffari   PetscScalar *x;
411ca3ac6ddSLeila Ghaffari   PetscErrorCode ierr;
412ca3ac6ddSLeila Ghaffari 
413ca3ac6ddSLeila Ghaffari   PetscFunctionBeginUser;
414ca3ac6ddSLeila Ghaffari   // Composite Operaters
415ca3ac6ddSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
4161e150236SLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply, op_applyVol); // Apply a Sub-Operator for the volume
417ca3ac6ddSLeila Ghaffari 
418*7659d40cSLeila Ghaffari   if (wind_type == ADVECTION_WIND_TRANSLATION) {
419*7659d40cSLeila Ghaffari     bc->nwall = 0;
420*7659d40cSLeila Ghaffari     bc->nslip[0] = bc->nslip[1] = bc->nslip[2] = 0;
4211e150236SLeila Ghaffari     ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
4221e150236SLeila Ghaffari     ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr);
4231e150236SLeila Ghaffari     ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr);
4241e150236SLeila Ghaffari     ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr);
4251e150236SLeila Ghaffari     CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x);
426ca3ac6ddSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
427*7659d40cSLeila Ghaffari     ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
428*7659d40cSLeila Ghaffari     if (dim == 2) nFace = 4;
429*7659d40cSLeila Ghaffari     if (dim == 3) nFace = 6;
4309fe13df9SLeila Ghaffari 
431*7659d40cSLeila Ghaffari     // Create CEED Operator for each boundary face
432*7659d40cSLeila Ghaffari     for (CeedInt i=0; i<nFace; i++) {
433*7659d40cSLeila Ghaffari       ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, i+1, numP_Sur,
434*7659d40cSLeila Ghaffari                                      numQ_Sur, qdatasizeSur, &restrictqSur[i], &restrictxSur[i],
435*7659d40cSLeila Ghaffari                                      &restrictqdiSur[i]); CHKERRQ(ierr);
436*7659d40cSLeila Ghaffari       // Create the CEED vectors that will be needed in Boundary setup
437*7659d40cSLeila Ghaffari       CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
438*7659d40cSLeila Ghaffari       CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur[i]);
439*7659d40cSLeila Ghaffari       // Create the operator that builds the quadrature data for the Boundary operator
440*7659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
441*7659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur, CEED_VECTOR_ACTIVE);
442*7659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
443ca3ac6ddSLeila Ghaffari                            basisxSur, CEED_VECTOR_NONE);
444*7659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
445ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
446*7659d40cSLeila Ghaffari       // Create Boundary operator
447*7659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]);
448*7659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
449*7659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i],
450*7659d40cSLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
451*7659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur, xcorners);
452*7659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur, CEED_VECTOR_ACTIVE);
453*7659d40cSLeila Ghaffari       // Apply CEED operator for Boundary setup
454*7659d40cSLeila Ghaffari       CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i], CEED_REQUEST_IMMEDIATE);
455*7659d40cSLeila Ghaffari       // Apply Sub-Operator for the Boundary
456*7659d40cSLeila Ghaffari       CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]);
4579fe13df9SLeila Ghaffari     }
4581e150236SLeila Ghaffari     CeedVectorDestroy(&xcorners);
459ca3ac6ddSLeila Ghaffari   }
460ca3ac6ddSLeila Ghaffari   PetscFunctionReturn(0);
461ca3ac6ddSLeila Ghaffari }
462ca3ac6ddSLeila Ghaffari 
463ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
464ccaff030SJeremy L Thompson   PetscErrorCode ierr;
465ccaff030SJeremy L Thompson   PetscInt m;
466ccaff030SJeremy L Thompson 
467ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
468ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
469ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
470ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
471ccaff030SJeremy L Thompson }
472ccaff030SJeremy L Thompson 
473ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
474ccaff030SJeremy L Thompson   PetscErrorCode ierr;
475ccaff030SJeremy L Thompson   PetscInt mceed, mpetsc;
476ccaff030SJeremy L Thompson   PetscScalar *a;
477ccaff030SJeremy L Thompson 
478ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
479ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
480ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
481ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
48284d34d69SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
48384d34d69SLeila Ghaffari                                   mpetsc, mceed);
484ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
485ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
486ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
487ccaff030SJeremy L Thompson }
488ccaff030SJeremy L Thompson 
489ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
490ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
491ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
492ccaff030SJeremy L Thompson   PetscErrorCode ierr;
493ccaff030SJeremy L Thompson   Vec Qbc;
494ccaff030SJeremy L Thompson 
495ccaff030SJeremy L Thompson   PetscFunctionBegin;
496ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
497ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
499ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
500ccaff030SJeremy L Thompson }
501ccaff030SJeremy L Thompson 
502ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
503ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
504ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
505ccaff030SJeremy L Thompson   PetscErrorCode ierr;
506ccaff030SJeremy L Thompson   User user = *(User *)userData;
507ccaff030SJeremy L Thompson   PetscScalar *q, *g;
508ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
509ccaff030SJeremy L Thompson 
510ccaff030SJeremy L Thompson   // Global-to-local
511ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
512ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
513ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
514ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
515ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
516ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
517ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
518ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
519ccaff030SJeremy L Thompson 
520ccaff030SJeremy L Thompson   // Ceed Vectors
521ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
522ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
523ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
524ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
525ccaff030SJeremy L Thompson 
526ccaff030SJeremy L Thompson   // Apply CEED operator
527ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
528ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
529ccaff030SJeremy L Thompson 
530ccaff030SJeremy L Thompson   // Restore vectors
531ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
532ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
533ccaff030SJeremy L Thompson 
534ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
535ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
536ccaff030SJeremy L Thompson 
537ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
538ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
539ccaff030SJeremy L Thompson   CHKERRQ(ierr);
540ccaff030SJeremy L Thompson 
541ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
542ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
543ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
544ccaff030SJeremy L Thompson }
545ccaff030SJeremy L Thompson 
546ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
547ccaff030SJeremy L Thompson                                    void *userData) {
548ccaff030SJeremy L Thompson   PetscErrorCode ierr;
549ccaff030SJeremy L Thompson   User user = *(User *)userData;
550ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
551ccaff030SJeremy L Thompson   PetscScalar *g;
552ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
553ccaff030SJeremy L Thompson 
554ccaff030SJeremy L Thompson   // Global-to-local
555ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
556ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
557ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
558ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
559ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
560ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
561ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
562ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
563ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
564ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
565ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
566ccaff030SJeremy L Thompson 
567ccaff030SJeremy L Thompson   // Ceed Vectors
568ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
569ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
570ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
571ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
572ccaff030SJeremy L Thompson                      (PetscScalar *)q);
573ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
574ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
575ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
576ccaff030SJeremy L Thompson 
577ccaff030SJeremy L Thompson   // Apply CEED operator
578ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
579ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
580ccaff030SJeremy L Thompson 
581ccaff030SJeremy L Thompson   // Restore vectors
582ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
583ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
584ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
585ccaff030SJeremy L Thompson 
586ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
587ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
588ccaff030SJeremy L Thompson 
589ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
590ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
591ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
592ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
593ccaff030SJeremy L Thompson }
594ccaff030SJeremy L Thompson 
595ccaff030SJeremy L Thompson // User provided TS Monitor
596ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
597ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
598ccaff030SJeremy L Thompson   User user = ctx;
599ccaff030SJeremy L Thompson   Vec Qloc;
600ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
601ccaff030SJeremy L Thompson   PetscViewer viewer;
602ccaff030SJeremy L Thompson   PetscErrorCode ierr;
603ccaff030SJeremy L Thompson 
604ccaff030SJeremy L Thompson   // Set up output
605ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
606ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
607ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
608ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
609ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
610ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
611ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
612ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
613ccaff030SJeremy L Thompson 
614ccaff030SJeremy L Thompson   // Output
615ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
616ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
617ccaff030SJeremy L Thompson   CHKERRQ(ierr);
618ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
619ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
620ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
6219d801c56SJed Brown   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
622ccaff030SJeremy L Thompson   if (user->dmviz) {
623ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
624ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
625ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
626ccaff030SJeremy L Thompson 
627ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
628ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
629ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
630ccaff030SJeremy L Thompson     CHKERRQ(ierr);
631ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
632ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
633ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
634ccaff030SJeremy L Thompson     CHKERRQ(ierr);
635ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
636ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
637ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
638ccaff030SJeremy L Thompson     CHKERRQ(ierr);
639ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
640ccaff030SJeremy L Thompson                               filepath_refined,
641ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
642ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
643ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
644ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
645ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
646ccaff030SJeremy L Thompson   }
647ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
648ccaff030SJeremy L Thompson 
649ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
650ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
651ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
652ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
653ccaff030SJeremy L Thompson   CHKERRQ(ierr);
654ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
655ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
656ccaff030SJeremy L Thompson 
657ccaff030SJeremy L Thompson   // Save time stamp
658ccaff030SJeremy L Thompson   // Dimensionalize time back
659ccaff030SJeremy L Thompson   time /= user->units->second;
660ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
661ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
662ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
663ccaff030SJeremy L Thompson   CHKERRQ(ierr);
664ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
665ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
666ccaff030SJeremy L Thompson   #else
667ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
668ccaff030SJeremy L Thompson   #endif
669ccaff030SJeremy L Thompson   CHKERRQ(ierr);
670ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
671ccaff030SJeremy L Thompson 
672ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
673ccaff030SJeremy L Thompson }
674ccaff030SJeremy L Thompson 
67584d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
676ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
677ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
678ccaff030SJeremy L Thompson   PetscErrorCode ierr;
679ccaff030SJeremy L Thompson   CeedVector multlvec;
680ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
681ccaff030SJeremy L Thompson 
682ccaff030SJeremy L Thompson   ctxSetup->time = time;
683ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
684ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
685ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
686ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
687ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
688ccaff030SJeremy L Thompson 
689ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
690ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
691ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
692ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
693ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
694ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
695ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
696ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
697ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
698ccaff030SJeremy L Thompson   CHKERRQ(ierr);
699ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
700ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
701ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
702ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
703ccaff030SJeremy L Thompson 
704ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
705ccaff030SJeremy L Thompson }
706ccaff030SJeremy L Thompson 
707ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
708ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
709ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
710ccaff030SJeremy L Thompson   PetscErrorCode ierr;
711ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
712ccaff030SJeremy L Thompson   CeedOperator op_mass;
713ccaff030SJeremy L Thompson   CeedVector mceed;
714ccaff030SJeremy L Thompson   Vec Mloc;
715ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
716ccaff030SJeremy L Thompson 
717ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
718ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
719ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
720ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
721ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
722ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
723ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
724ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
725ccaff030SJeremy L Thompson 
726ccaff030SJeremy L Thompson   // Create the mass operator
727ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
728ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
729ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
730ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
731ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
732ccaff030SJeremy L Thompson 
733ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
734ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
735ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
736ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
737ccaff030SJeremy L Thompson 
738ccaff030SJeremy L Thompson   {
739ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
740ccaff030SJeremy L Thompson     CeedVector onesvec;
741ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
742ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
743ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
744ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
745ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
746ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
747ccaff030SJeremy L Thompson   }
748ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
749ccaff030SJeremy L Thompson 
750ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
751ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
752ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
753ccaff030SJeremy L Thompson 
754ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
755ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
756ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
757ccaff030SJeremy L Thompson }
758ccaff030SJeremy L Thompson 
75984d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
760ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
761ccaff030SJeremy L Thompson   PetscErrorCode ierr;
762ccaff030SJeremy L Thompson 
763ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
764ccaff030SJeremy L Thompson   {
765ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
766ccaff030SJeremy L Thompson     PetscFE fe;
767ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
768ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
769ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
77032ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
771ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
772ccaff030SJeremy L Thompson     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
773ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
77407af6069Svaleriabarra     {
77507af6069Svaleriabarra       PetscInt comps[1] = {1};
77607af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
77707af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
77807af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
77907af6069Svaleriabarra       comps[0] = 2;
78007af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
78107af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
78207af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
78307af6069Svaleriabarra       comps[0] = 3;
78407af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
78507af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
78607af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
78707af6069Svaleriabarra     }
78884d34d69SLeila Ghaffari     if (bc->userbc == PETSC_TRUE) {
78984d34d69SLeila Ghaffari       for (PetscInt c = 0; c < 3; c++) {
79084d34d69SLeila Ghaffari         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
79184d34d69SLeila Ghaffari           for (PetscInt w = 0; w < bc->nwall; w++) {
79284d34d69SLeila Ghaffari             if (bc->slips[c][s] == bc->walls[w])
79384d34d69SLeila Ghaffari               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
79484d34d69SLeila Ghaffari                        "Boundary condition already set on face %D!\n", bc->walls[w]);
79584d34d69SLeila Ghaffari 
79684d34d69SLeila Ghaffari           }
79784d34d69SLeila Ghaffari         }
79884d34d69SLeila Ghaffari       }
79984d34d69SLeila Ghaffari     }
80084d34d69SLeila Ghaffari     // Wall boundary conditions are zero energy density and zero flux for
80184d34d69SLeila Ghaffari     //   velocity in advection/advection2d, and zero velocity and zero flux
80284d34d69SLeila Ghaffari     //   for mass density and energy density in density_current
80384d34d69SLeila Ghaffari     {
80484d34d69SLeila Ghaffari       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
80584d34d69SLeila Ghaffari         PetscInt comps[1] = {4};
80684d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
80784d34d69SLeila Ghaffari                              1, comps, (void(*)(void))problem->bc,
80884d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
80984d34d69SLeila Ghaffari       } else if (problem->bc == Exact_DC) {
81084d34d69SLeila Ghaffari         PetscInt comps[3] = {1, 2, 3};
81184d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
81284d34d69SLeila Ghaffari                              3, comps, (void(*)(void))problem->bc,
81384d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
81484d34d69SLeila Ghaffari       } else
81584d34d69SLeila Ghaffari         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
81684d34d69SLeila Ghaffari                 "Undefined boundary conditions for this problem");
81784d34d69SLeila Ghaffari     }
818ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
819ccaff030SJeremy L Thompson     CHKERRQ(ierr);
820ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
821ccaff030SJeremy L Thompson   }
822ccaff030SJeremy L Thompson   {
823ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
824ccaff030SJeremy L Thompson     PetscSection section;
825ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
826ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
827ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
828ccaff030SJeremy L Thompson     CHKERRQ(ierr);
829ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
830ccaff030SJeremy L Thompson     CHKERRQ(ierr);
831ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
832ccaff030SJeremy L Thompson     CHKERRQ(ierr);
833ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
834ccaff030SJeremy L Thompson     CHKERRQ(ierr);
835ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
836ccaff030SJeremy L Thompson     CHKERRQ(ierr);
837ccaff030SJeremy L Thompson   }
838ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
839ccaff030SJeremy L Thompson }
840ccaff030SJeremy L Thompson 
841ccaff030SJeremy L Thompson int main(int argc, char **argv) {
842ccaff030SJeremy L Thompson   PetscInt ierr;
843ccaff030SJeremy L Thompson   MPI_Comm comm;
84484d34d69SLeila Ghaffari   DM dm, dmcoord, dmviz;
845ccaff030SJeremy L Thompson   Mat interpviz;
846ccaff030SJeremy L Thompson   TS ts;
847ccaff030SJeremy L Thompson   TSAdapt adapt;
848ccaff030SJeremy L Thompson   User user;
849ccaff030SJeremy L Thompson   Units units;
850ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
85184d34d69SLeila Ghaffari   PetscInt localNelemVol, lnodes, gnodes, steps;
852ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
853ccaff030SJeremy L Thompson   PetscMPIInt rank;
854ccaff030SJeremy L Thompson   PetscScalar ftime;
855ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
856ccaff030SJeremy L Thompson   Ceed ceed;
85784d34d69SLeila Ghaffari   CeedInt numP, numQ;
858cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
85984d34d69SLeila Ghaffari   CeedBasis basisx, basisxc, basisq;
86084d34d69SLeila Ghaffari   CeedElemRestriction restrictx, restrictq, restrictqdi;
861cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
862cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
863ccaff030SJeremy L Thompson   CeedScalar Rd;
86484d34d69SLeila Ghaffari   CeedMemType memtyperequested;
865ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
866ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
867ccaff030SJeremy L Thompson   problemType problemChoice;
868ccaff030SJeremy L Thompson   problemData *problem = NULL;
8691184866aSLeila Ghaffari   WindType wind_type;
870ccaff030SJeremy L Thompson   StabilizationType stab;
87184d34d69SLeila Ghaffari   testType testChoice;
87284d34d69SLeila Ghaffari   testData *test = NULL;
87384d34d69SLeila Ghaffari   PetscBool implicit;
874cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
875ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
87684d34d69SLeila Ghaffari     .nslip = {2, 2, 2},
87784d34d69SLeila Ghaffari     .slips = {{5, 6}, {3, 4}, {1, 2}}
878ccaff030SJeremy L Thompson   };
879ccaff030SJeremy L Thompson   double start, cpu_time_used;
88084d34d69SLeila Ghaffari   // Check PETSc CUDA support
88184d34d69SLeila Ghaffari   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
88284d34d69SLeila Ghaffari   // *INDENT-OFF*
88384d34d69SLeila Ghaffari   #ifdef PETSC_HAVE_CUDA
88484d34d69SLeila Ghaffari   petschavecuda = PETSC_TRUE;
88584d34d69SLeila Ghaffari   #else
88684d34d69SLeila Ghaffari   petschavecuda = PETSC_FALSE;
88784d34d69SLeila Ghaffari   #endif
88884d34d69SLeila Ghaffari   // *INDENT-ON*
889ccaff030SJeremy L Thompson 
890ccaff030SJeremy L Thompson   // Create the libCEED contexts
891ccaff030SJeremy L Thompson   PetscScalar meter      = 1e-2;     // 1 meter in scaled length units
892ccaff030SJeremy L Thompson   PetscScalar second     = 1e-2;     // 1 second in scaled time units
893ccaff030SJeremy L Thompson   PetscScalar kilogram   = 1e-6;     // 1 kilogram in scaled mass units
894ccaff030SJeremy L Thompson   PetscScalar Kelvin     = 1;        // 1 Kelvin in scaled temperature units
895ccaff030SJeremy L Thompson   CeedScalar theta0      = 300.;     // K
896ccaff030SJeremy L Thompson   CeedScalar thetaC      = -15.;     // K
897ccaff030SJeremy L Thompson   CeedScalar P0          = 1.e5;     // Pa
89881f92cf0SLeila Ghaffari   CeedScalar P_wind      = 1.5e5;    // Pa
89981f92cf0SLeila Ghaffari   CeedScalar rho_wind    = 1.2;      // Kg/m^3
900ccaff030SJeremy L Thompson   CeedScalar N           = 0.01;     // 1/s
901ccaff030SJeremy L Thompson   CeedScalar cv          = 717.;     // J/(kg K)
902ccaff030SJeremy L Thompson   CeedScalar cp          = 1004.;    // J/(kg K)
903ccaff030SJeremy L Thompson   CeedScalar g           = 9.81;     // m/s^2
904ccaff030SJeremy L Thompson   CeedScalar lambda      = -2./3.;   // -
905ccaff030SJeremy L Thompson   CeedScalar mu          = 75.;      // Pa s, dynamic viscosity
906ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
907ccaff030SJeremy L Thompson   CeedScalar k           = 0.02638;  // W/(m K)
908ccaff030SJeremy L Thompson   CeedScalar CtauS       = 0.;       // dimensionless
909ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;       // [0,1]
910ccaff030SJeremy L Thompson   PetscScalar lx         = 8000.;    // m
911ccaff030SJeremy L Thompson   PetscScalar ly         = 8000.;    // m
912ccaff030SJeremy L Thompson   PetscScalar lz         = 4000.;    // m
913ccaff030SJeremy L Thompson   CeedScalar rc          = 1000.;    // m (Radius of bubble)
914ccaff030SJeremy L Thompson   PetscScalar resx       = 1000.;    // m (resolution in x)
915ccaff030SJeremy L Thompson   PetscScalar resy       = 1000.;    // m (resolution in y)
916ccaff030SJeremy L Thompson   PetscScalar resz       = 1000.;    // m (resolution in z)
917ccaff030SJeremy L Thompson   PetscInt outputfreq    = 10;       // -
918ccaff030SJeremy L Thompson   PetscInt contsteps     = 0;        // -
91984d34d69SLeila Ghaffari   PetscInt degree        = 1;        // -
92084d34d69SLeila Ghaffari   PetscInt qextra        = 2;        // -
921ea6e0f84SLeila Ghaffari   PetscInt qextraSur     = 2;        // -
9221184866aSLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0}, wind[3] = {1., 0, 0};
923ccaff030SJeremy L Thompson 
924ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
925ccaff030SJeremy L Thompson   if (ierr) return ierr;
926ccaff030SJeremy L Thompson 
927ccaff030SJeremy L Thompson   // Allocate PETSc context
928ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
929ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
930ccaff030SJeremy L Thompson 
931ccaff030SJeremy L Thompson   // Parse command line options
932ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
933ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
934ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
935ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
936ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
937ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
93884d34d69SLeila Ghaffari   testChoice = TEST_NONE;
93984d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-test", "Run tests", NULL,
94084d34d69SLeila Ghaffari                           testTypes, (PetscEnum)testChoice,
94184d34d69SLeila Ghaffari                           (PetscEnum *)&testChoice,
94284d34d69SLeila Ghaffari                           NULL); CHKERRQ(ierr);
94384d34d69SLeila Ghaffari   test = &testOptions[testChoice];
944ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
945ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
946ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
947ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
948ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
9491184866aSLeila Ghaffari   ierr = PetscOptionsEnum("-problem_advection_wind", "Wind type in Advection",
9501184866aSLeila Ghaffari                           NULL, WindTypes, (PetscEnum)(wind_type = ADVECTION_WIND_ROTATION),
9511184866aSLeila Ghaffari                           (PetscEnum *)&wind_type, NULL); CHKERRQ(ierr);
9521184866aSLeila Ghaffari   PetscInt n = problem->dim;
95382c09b01SLeila Ghaffari   PetscBool userWind;
9541184866aSLeila Ghaffari   ierr = PetscOptionsRealArray("-problem_advection_wind_translation", "Constant wind vector",
95582c09b01SLeila Ghaffari                                NULL, wind, &n, &userWind); CHKERRQ(ierr);
95682c09b01SLeila Ghaffari   if (wind_type == ADVECTION_WIND_ROTATION && userWind) {
95782c09b01SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -problem_advection_wind_translation only with -problem_advection_wind translation\n");
95882c09b01SLeila Ghaffari     CHKERRQ(ierr);
9591184866aSLeila Ghaffari   }
960ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
961ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
962ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
963ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
964ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
965ccaff030SJeremy L Thompson   CHKERRQ(ierr);
96684d34d69SLeila Ghaffari   if (!implicit && stab != STAB_NONE) {
96784d34d69SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
96884d34d69SLeila Ghaffari     CHKERRQ(ierr);
96984d34d69SLeila Ghaffari   }
970ccaff030SJeremy L Thompson   {
9717573aee6SJed Brown     PetscInt len;
9727573aee6SJed Brown     PetscBool flg;
973ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
974ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
975ccaff030SJeremy L Thompson                                 NULL, bc.walls,
976ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
977ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
9787573aee6SJed Brown     if (flg) {
9797573aee6SJed Brown       bc.nwall = len;
9807573aee6SJed Brown       // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
9817573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
9827573aee6SJed Brown     }
983ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
984ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
985ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
986ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
987ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
988ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
989ccaff030SJeremy L Thompson                                    &len), &flg);
990ccaff030SJeremy L Thompson       CHKERRQ(ierr);
99184d34d69SLeila Ghaffari       if (flg) {
99284d34d69SLeila Ghaffari         bc.nslip[j] = len;
99384d34d69SLeila Ghaffari         bc.userbc = PETSC_TRUE;
99484d34d69SLeila Ghaffari       }
995ccaff030SJeremy L Thompson     }
996ccaff030SJeremy L Thompson   }
997cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
998cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
999cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
1000ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1001ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
1002ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
1003ccaff030SJeremy L Thompson   meter = fabs(meter);
1004ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
1005ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
1006ccaff030SJeremy L Thompson   second = fabs(second);
1007ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
1008ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
1009ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
1010ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
1011ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
1012ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
1013ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
1014ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
1015ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
1016ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
1017ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
1018ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
1019ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
102081f92cf0SLeila Ghaffari   ierr = PetscOptionsScalar("-P_wind", "Inflow wind pressure",
102181f92cf0SLeila Ghaffari                             NULL, P_wind, &P_wind, NULL); CHKERRQ(ierr);
102281f92cf0SLeila Ghaffari   ierr = PetscOptionsScalar("-rho_wind", "Inflow wind density",
102381f92cf0SLeila Ghaffari                             NULL, rho_wind, &rho_wind, NULL); CHKERRQ(ierr);
1024ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
1025ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
1026ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
1027ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
1028ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
1029ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
1030ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
1031ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
1032ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
1033ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
1034ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1035ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1036ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1037ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1038ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
1039ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
1040ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
1041ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
104284d34d69SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
104384d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
104484d34d69SLeila Ghaffari                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
104584d34d69SLeila Ghaffari     CHKERRQ(ierr);
104684d34d69SLeila Ghaffari   }
1047ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
1048ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
1049ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
1050ccaff030SJeremy L Thompson   CHKERRQ(ierr);
105184d34d69SLeila Ghaffari   if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) {
105284d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
105384d34d69SLeila Ghaffari                        "Warning! Problem density_current does not support -CtauS or -strong_form\n");
105484d34d69SLeila Ghaffari     CHKERRQ(ierr);
105584d34d69SLeila Ghaffari   }
1056ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
1057ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
1058ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
1059ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
1060ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
1061ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
1062ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
1063ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
1064ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
1065ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
1066ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
1067ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
1068ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
1069ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
107082c09b01SLeila Ghaffari   n = problem->dim;
1071ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
1072ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
1073ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
1074ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
1075ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
1076ccaff030SJeremy L Thompson   n = problem->dim;
1077ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
1078ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
1079ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
1080ccaff030SJeremy L Thompson   {
1081ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
1082ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
1083ccaff030SJeremy L Thompson     if (norm > 0) {
1084ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
1085ccaff030SJeremy L Thompson     }
1086ccaff030SJeremy L Thompson   }
1087ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
1088ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
1089ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
1090ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
1091ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
109284d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
109384d34d69SLeila Ghaffari                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
109484d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
109584d34d69SLeila Ghaffari                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
109681f92cf0SLeila Ghaffari   PetscBool userQextraSur;
109781f92cf0SLeila Ghaffari   ierr = PetscOptionsInt("-qextra_boundary", "Number of extra quadrature points on in/outflow faces",
109881f92cf0SLeila Ghaffari                          NULL, qextraSur, &qextraSur, &userQextraSur); CHKERRQ(ierr);
109981f92cf0SLeila Ghaffari   if ((wind_type == ADVECTION_WIND_ROTATION || problemChoice == NS_DENSITY_CURRENT) && userQextraSur) {
110081f92cf0SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -qextra_boundary only with -problem_advection_wind translation\n");
110181f92cf0SLeila Ghaffari     CHKERRQ(ierr);
110281f92cf0SLeila Ghaffari   }
110384d34d69SLeila Ghaffari   ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr);
1104ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
1105ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
1106ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
110784d34d69SLeila Ghaffari   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
110884d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-memtype",
110984d34d69SLeila Ghaffari                           "CEED MemType requested", NULL,
111084d34d69SLeila Ghaffari                           memTypes, (PetscEnum)memtyperequested,
111184d34d69SLeila Ghaffari                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
111284d34d69SLeila Ghaffari   CHKERRQ(ierr);
1113ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1114ccaff030SJeremy L Thompson 
1115ccaff030SJeremy L Thompson   // Define derived units
1116ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
1117ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1118ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
1119ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
1120ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
1121ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
1122ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
1123ccaff030SJeremy L Thompson 
1124ccaff030SJeremy L Thompson   // Scale variables to desired units
1125ccaff030SJeremy L Thompson   theta0 *= Kelvin;
1126ccaff030SJeremy L Thompson   thetaC *= Kelvin;
1127ccaff030SJeremy L Thompson   P0 *= Pascal;
112881f92cf0SLeila Ghaffari   P_wind *= Pascal;
112981f92cf0SLeila Ghaffari   rho_wind *= kgpercubicm;
1130ccaff030SJeremy L Thompson   N *= (1./second);
1131ccaff030SJeremy L Thompson   cv *= JperkgK;
1132ccaff030SJeremy L Thompson   cp *= JperkgK;
1133ccaff030SJeremy L Thompson   Rd = cp - cv;
1134ccaff030SJeremy L Thompson   g *= mpersquareds;
1135ccaff030SJeremy L Thompson   mu *= Pascal * second;
1136ccaff030SJeremy L Thompson   k *= WpermK;
1137ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
1138ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
1139ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
1140ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
1141ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
1142ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
1143ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
1144ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
1145ccaff030SJeremy L Thompson 
1146ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
1147cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
1148ccaff030SJeremy L Thompson   // Set up the libCEED context
1149ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
1150ccaff030SJeremy L Thompson     .theta0 = theta0,
1151ccaff030SJeremy L Thompson     .thetaC = thetaC,
1152ccaff030SJeremy L Thompson     .P0 = P0,
1153ccaff030SJeremy L Thompson     .N = N,
1154ccaff030SJeremy L Thompson     .cv = cv,
1155ccaff030SJeremy L Thompson     .cp = cp,
1156ccaff030SJeremy L Thompson     .Rd = Rd,
1157ccaff030SJeremy L Thompson     .g = g,
1158ccaff030SJeremy L Thompson     .rc = rc,
1159ccaff030SJeremy L Thompson     .lx = lx,
1160ccaff030SJeremy L Thompson     .ly = ly,
1161ccaff030SJeremy L Thompson     .lz = lz,
1162ccaff030SJeremy L Thompson     .center[0] = center[0],
1163ccaff030SJeremy L Thompson     .center[1] = center[1],
1164ccaff030SJeremy L Thompson     .center[2] = center[2],
1165ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
1166ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
1167ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
11681184866aSLeila Ghaffari     .wind[0] = wind[0],
11691184866aSLeila Ghaffari     .wind[1] = wind[1],
11701184866aSLeila Ghaffari     .wind[2] = wind[2],
1171ccaff030SJeremy L Thompson     .time = 0,
11721184866aSLeila Ghaffari     .wind_type = wind_type,
1173ccaff030SJeremy L Thompson   };
1174ccaff030SJeremy L Thompson 
117584d34d69SLeila Ghaffari   // Create the mesh
1176ccaff030SJeremy L Thompson   {
1177ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
1178ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
117984d34d69SLeila Ghaffari                                NULL, PETSC_TRUE, &dm);
1180ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1181ccaff030SJeremy L Thompson   }
118284d34d69SLeila Ghaffari 
118384d34d69SLeila Ghaffari   // Distribute the mesh over processes
118484d34d69SLeila Ghaffari   {
1185ccaff030SJeremy L Thompson     DM               dmDist = NULL;
1186ccaff030SJeremy L Thompson     PetscPartitioner part;
1187ccaff030SJeremy L Thompson 
1188ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1189ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1190ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1191ccaff030SJeremy L Thompson     if (dmDist) {
1192ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1193ccaff030SJeremy L Thompson       dm  = dmDist;
1194ccaff030SJeremy L Thompson     }
1195ccaff030SJeremy L Thompson   }
1196ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
1197ccaff030SJeremy L Thompson 
119884d34d69SLeila Ghaffari   // Setup DM
1199ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
1200ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
120184d34d69SLeila Ghaffari   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
120284d34d69SLeila Ghaffari 
120384d34d69SLeila Ghaffari   // Refine DM for high-order viz
1204ccaff030SJeremy L Thompson   dmviz = NULL;
1205ccaff030SJeremy L Thompson   interpviz = NULL;
1206ccaff030SJeremy L Thompson   if (viz_refine) {
1207ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
1208ff6701fcSJed Brown 
1209ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
1210ff6701fcSJed Brown     dmhierarchy[0] = dm;
121184d34d69SLeila Ghaffari     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
1212ff6701fcSJed Brown       Mat interp_next;
1213ff6701fcSJed Brown 
1214ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
1215ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1216ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
1217ff6701fcSJed Brown       d = (d + 1) / 2;
1218ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
1219ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
1220ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
1221ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
1222ff6701fcSJed Brown       if (!i) interpviz = interp_next;
1223ff6701fcSJed Brown       else {
1224ff6701fcSJed Brown         Mat C;
1225ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1226ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1227ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1228ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1229ff6701fcSJed Brown         interpviz = C;
1230ff6701fcSJed Brown       }
1231ff6701fcSJed Brown     }
1232cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1233ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1234cb3e2689Svaleriabarra     }
1235ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1236ccaff030SJeremy L Thompson   }
1237ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1238ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1239ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1240ccaff030SJeremy L Thompson   lnodes /= ncompq;
1241ccaff030SJeremy L Thompson 
124284d34d69SLeila Ghaffari   // Initialize CEED
124384d34d69SLeila Ghaffari   CeedInit(ceedresource, &ceed);
124484d34d69SLeila Ghaffari   // Set memtype
124584d34d69SLeila Ghaffari   CeedMemType memtypebackend;
124684d34d69SLeila Ghaffari   CeedGetPreferredMemType(ceed, &memtypebackend);
124784d34d69SLeila Ghaffari   // Check memtype compatibility
124884d34d69SLeila Ghaffari   if (!setmemtyperequest)
124984d34d69SLeila Ghaffari     memtyperequested = memtypebackend;
125084d34d69SLeila Ghaffari   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
125184d34d69SLeila Ghaffari     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
125284d34d69SLeila Ghaffari              "PETSc was not built with CUDA. "
125384d34d69SLeila Ghaffari              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
125484d34d69SLeila Ghaffari 
125584d34d69SLeila Ghaffari   // Set number of 1D nodes and quadrature points
125684d34d69SLeila Ghaffari   numP = degree + 1;
125784d34d69SLeila Ghaffari   numQ = numP + qextra;
125884d34d69SLeila Ghaffari 
125984d34d69SLeila Ghaffari     // Print summary
126084d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1261ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1262ccaff030SJeremy L Thompson     int comm_size;
1263ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1264ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1265ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
126684d34d69SLeila Ghaffari     gnodes = gdofs/ncompq;
1267ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1268ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1269ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
127084d34d69SLeila Ghaffari     const char *usedresource;
127184d34d69SLeila Ghaffari     CeedGetResource(ceed, &usedresource);
1272ccaff030SJeremy L Thompson 
127384d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
127484d34d69SLeila Ghaffari                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
127584d34d69SLeila Ghaffari                        "  rank(s)                              : %d\n"
127684d34d69SLeila Ghaffari                        "  Problem:\n"
127784d34d69SLeila Ghaffari                        "    Problem Name                       : %s\n"
127884d34d69SLeila Ghaffari                        "    Stabilization                      : %s\n"
127984d34d69SLeila Ghaffari                        "  PETSc:\n"
128084d34d69SLeila Ghaffari                        "    Box Faces                          : %s\n"
128184d34d69SLeila Ghaffari                        "  libCEED:\n"
128284d34d69SLeila Ghaffari                        "    libCEED Backend                    : %s\n"
128384d34d69SLeila Ghaffari                        "    libCEED Backend MemType            : %s\n"
128484d34d69SLeila Ghaffari                        "    libCEED User Requested MemType     : %s\n"
128584d34d69SLeila Ghaffari                        "  Mesh:\n"
128684d34d69SLeila Ghaffari                        "    Number of 1D Basis Nodes (P)       : %d\n"
128784d34d69SLeila Ghaffari                        "    Number of 1D Quadrature Points (Q) : %d\n"
128884d34d69SLeila Ghaffari                        "    Global DoFs                        : %D\n"
128984d34d69SLeila Ghaffari                        "    Owned DoFs                         : %D\n"
129084d34d69SLeila Ghaffari                        "    DoFs per node                      : %D\n"
129184d34d69SLeila Ghaffari                        "    Global nodes                       : %D\n"
129284d34d69SLeila Ghaffari                        "    Owned nodes                        : %D\n",
129384d34d69SLeila Ghaffari                        comm_size, problemTypes[problemChoice],
129484d34d69SLeila Ghaffari                        StabilizationTypes[stab], box_faces_str, usedresource,
129584d34d69SLeila Ghaffari                        CeedMemTypes[memtypebackend],
129684d34d69SLeila Ghaffari                        (setmemtyperequest) ?
129784d34d69SLeila Ghaffari                        CeedMemTypes[memtyperequested] : "none",
129884d34d69SLeila Ghaffari                        numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes);
129984d34d69SLeila Ghaffari     CHKERRQ(ierr);
13000c6c0b13SLeila Ghaffari   }
13010c6c0b13SLeila Ghaffari 
1302ccaff030SJeremy L Thompson   // Set up global mass vector
1303ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
1304ccaff030SJeremy L Thompson 
130584d34d69SLeila Ghaffari   // Set up libCEED
1306ccaff030SJeremy L Thompson   // CEED Bases
1307ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
130884d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
130984d34d69SLeila Ghaffari                                   &basisq);
131084d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
131184d34d69SLeila Ghaffari                                   &basisx);
131284d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
131384d34d69SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxc);
1314ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1315ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1316ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1317ccaff030SJeremy L Thompson 
1318ccaff030SJeremy L Thompson   // CEED Restrictions
13191e150236SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ,
132084d34d69SLeila Ghaffari                                  qdatasizeVol, &restrictq, &restrictx,
132184d34d69SLeila Ghaffari                                  &restrictqdi); CHKERRQ(ierr);
1322ccaff030SJeremy L Thompson 
1323ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1324ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1325ccaff030SJeremy L Thompson 
1326ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1327bd910870SLeila Ghaffari   CeedInt NqptsVol;
132884d34d69SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol);
132984d34d69SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictq, &localNelemVol);
13308b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
133184d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1332ccaff030SJeremy L Thompson 
1333ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1334ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1335ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1336ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1337ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
13388b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1339ccaff030SJeremy L Thompson 
1340ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1341ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1342ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1343ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1344ccaff030SJeremy L Thompson 
1345ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1346ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1347ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1348ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1349ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1350ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
13518b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1352ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1353ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1354ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1355ccaff030SJeremy L Thompson   }
1356ccaff030SJeremy L Thompson 
1357ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1358ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1359ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1360ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1361ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1362ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1363ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
13648b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1365ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1366ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1367ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1368ccaff030SJeremy L Thompson   }
1369ccaff030SJeremy L Thompson 
1370ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1371ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
137284d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1373ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
137484d34d69SLeila Ghaffari                        basisx, CEED_VECTOR_NONE);
137584d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdi,
1376ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1377ccaff030SJeremy L Thompson 
1378ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1379ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
138084d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
138184d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictq,
1382ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1383ccaff030SJeremy L Thompson 
138484d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
138584d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
138684d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1387ccaff030SJeremy L Thompson 
1388ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1389ccaff030SJeremy L Thompson     CeedOperator op;
1390ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
139184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
139284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
139384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
13948b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
139584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
139684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
139784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1398d3630711SJed Brown     user->op_rhs_vol = op;
1399ccaff030SJeremy L Thompson   }
1400ccaff030SJeremy L Thompson 
1401ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1402ccaff030SJeremy L Thompson     CeedOperator op;
1403ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
140484d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
140584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
140684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
140784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
14088b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
140984d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
141084d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
141184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1412d3630711SJed Brown     user->op_ifunction_vol = op;
1413ccaff030SJeremy L Thompson   }
1414ccaff030SJeremy L Thompson 
14156a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
14166a0edaf9SLeila Ghaffari   CeedInt height = 1;
14176a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
14181e150236SLeila Ghaffari   CeedInt numP_Sur = degree + 1;
14191e150236SLeila Ghaffari   CeedInt numQ_Sur = numP_Sur + qextraSur;
1420cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1421cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1422cfa64770SLeila Ghaffari   CeedInt NqptsSur;
1423*7659d40cSLeila Ghaffari   CeedQFunction qf_setupSur, qf_Sur;
1424cfa64770SLeila Ghaffari 
1425cfa64770SLeila Ghaffari   // CEED bases for the boundaries
14266a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
14276a0edaf9SLeila Ghaffari                                   &basisqSur);
14286a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
14296a0edaf9SLeila Ghaffari                                   &basisxSur);
14306a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
14316a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
14326a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
14336a0edaf9SLeila Ghaffari 
1434cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
14356a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
14366a0edaf9SLeila Ghaffari                               &qf_setupSur);
14376a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
14386a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
14396a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14406a0edaf9SLeila Ghaffari 
1441*7659d40cSLeila Ghaffari   // Creat Q-Function for Boundaries
1442*7659d40cSLeila Ghaffari   qf_Sur = NULL;
1443*7659d40cSLeila Ghaffari   if (problem->applySur) {
1444*7659d40cSLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur,
1445*7659d40cSLeila Ghaffari                                 problem->applySur_loc, &qf_Sur);
1446*7659d40cSLeila Ghaffari     CeedQFunctionAddInput(qf_Sur, "q", ncompq, CEED_EVAL_INTERP);
1447*7659d40cSLeila Ghaffari     CeedQFunctionAddInput(qf_Sur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
1448*7659d40cSLeila Ghaffari     CeedQFunctionAddInput(qf_Sur, "x", ncompx, CEED_EVAL_INTERP);
1449*7659d40cSLeila Ghaffari     CeedQFunctionAddOutput(qf_Sur, "v", ncompq, CEED_EVAL_INTERP);
14509fe13df9SLeila Ghaffari   }
14519fe13df9SLeila Ghaffari 
14529fe13df9SLeila Ghaffari   // Create CEED Operator for the whole domain
14539fe13df9SLeila Ghaffari   if (!implicit)
1454*7659d40cSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_rhs_vol, qf_Sur, qf_setupSur,
1455*7659d40cSLeila Ghaffari                                   height, numP_Sur, numQ_Sur, qdatasizeSur, NqptsSur, basisxSur,
1456*7659d40cSLeila Ghaffari                                   basisqSur, &user->op_rhs); CHKERRQ(ierr);
14579fe13df9SLeila Ghaffari   if (implicit)
1458*7659d40cSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_ifunction_vol, qf_Sur, qf_setupSur,
1459*7659d40cSLeila Ghaffari                                   height, numP_Sur, numQ_Sur, qdatasizeSur, NqptsSur, basisxSur,
1460*7659d40cSLeila Ghaffari                                   basisqSur, &user->op_ifunction); CHKERRQ(ierr);
1461*7659d40cSLeila Ghaffari   // Set up contex for QFunctions
1462ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1463ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1464ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1465ccaff030SJeremy L Thompson     .CtauS = CtauS,
1466ccaff030SJeremy L Thompson     .strong_form = strong_form,
1467ccaff030SJeremy L Thompson     .stabilization = stab,
1468ccaff030SJeremy L Thompson   };
1469*7659d40cSLeila Ghaffari   struct SurfaceContext_ ctxSurface = {
1470*7659d40cSLeila Ghaffari     .cv = cv,
1471*7659d40cSLeila Ghaffari     .cp = cp,
1472*7659d40cSLeila Ghaffari     .Rd = Rd,
1473*7659d40cSLeila Ghaffari     .P_wind = P_wind,
1474*7659d40cSLeila Ghaffari     .rho_wind = rho_wind,
1475*7659d40cSLeila Ghaffari     .strong_form = strong_form,
1476*7659d40cSLeila Ghaffari     .wind[0] = wind[0],
1477*7659d40cSLeila Ghaffari     .wind[1] = wind[1],
1478*7659d40cSLeila Ghaffari     .wind[2] = wind[2],
1479*7659d40cSLeila Ghaffari     .wind_type = wind_type,
1480*7659d40cSLeila Ghaffari     .implicit = implicit,
1481*7659d40cSLeila Ghaffari   };
1482ccaff030SJeremy L Thompson   switch (problemChoice) {
1483ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1484ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1485ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1486ccaff030SJeremy L Thompson           sizeof ctxNS);
1487ccaff030SJeremy L Thompson     break;
1488ccaff030SJeremy L Thompson   case NS_ADVECTION:
1489ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1490ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1491ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1492ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1493ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1494*7659d40cSLeila Ghaffari     if (qf_Sur) CeedQFunctionSetContext(qf_Sur, &ctxSurface, sizeof ctxSurface);
1495ccaff030SJeremy L Thompson   }
1496ccaff030SJeremy L Thompson 
1497ccaff030SJeremy L Thompson   // Set up PETSc context
1498ccaff030SJeremy L Thompson   // Set up units structure
1499ccaff030SJeremy L Thompson   units->meter = meter;
1500ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1501ccaff030SJeremy L Thompson   units->second = second;
1502ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1503ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1504ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1505ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1506ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1507ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1508ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1509ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1510ccaff030SJeremy L Thompson 
1511ccaff030SJeremy L Thompson   // Set up user structure
1512ccaff030SJeremy L Thompson   user->comm = comm;
1513ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1514ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1515ccaff030SJeremy L Thompson   user->units = units;
1516ccaff030SJeremy L Thompson   user->dm = dm;
1517ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1518ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1519ccaff030SJeremy L Thompson   user->ceed = ceed;
1520ccaff030SJeremy L Thompson 
15218b982baeSLeila Ghaffari   // Calculate qdata and ICs
1522ccaff030SJeremy L Thompson   // Set up state global and local vectors
1523ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1524ccaff030SJeremy L Thompson 
1525cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1526ccaff030SJeremy L Thompson 
1527ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1528ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
15298b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
153084d34d69SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1531ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1532ccaff030SJeremy L Thompson 
153384d34d69SLeila Ghaffari   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
153484d34d69SLeila Ghaffari                              &ctxSetup, 0.0); CHKERRQ(ierr);
1535ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1536ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1537ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1538ccaff030SJeremy L Thompson     // slower execution.
1539ccaff030SJeremy L Thompson     Vec Qbc;
1540ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1541ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1542ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1543ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1544ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1545ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1546ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
154784d34d69SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
154884d34d69SLeila Ghaffari     CHKERRQ(ierr);
1549ccaff030SJeremy L Thompson   }
1550ccaff030SJeremy L Thompson 
1551ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1552ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1553ccaff030SJeremy L Thompson   // Gather initial Q values
1554ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1555ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1556ccaff030SJeremy L Thompson     PetscViewer viewer;
1557ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1558ccaff030SJeremy L Thompson     // Read input
1559ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1560ccaff030SJeremy L Thompson                          user->outputfolder);
1561ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1562ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1563ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1564ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1565ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1566ccaff030SJeremy L Thompson   }
1567ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1568ccaff030SJeremy L Thompson 
1569ccaff030SJeremy L Thompson // Create and setup TS
1570ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1571ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1572ccaff030SJeremy L Thompson   if (implicit) {
1573ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1574ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1575ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1576ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1577ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1578ccaff030SJeremy L Thompson     }
1579ccaff030SJeremy L Thompson   } else {
1580ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL,
1581ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1582ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1583ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1584ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1585ccaff030SJeremy L Thompson   }
1586ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1587ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1588ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
158984d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1590ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1591ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1592ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1593ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1594ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1595ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
159684d34d69SLeila Ghaffari     if (testChoice == TEST_NONE) {
1597ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1598ccaff030SJeremy L Thompson     }
1599ccaff030SJeremy L Thompson   } else { // continue from time of last output
1600ccaff030SJeremy L Thompson     PetscReal time;
1601ccaff030SJeremy L Thompson     PetscInt count;
1602ccaff030SJeremy L Thompson     PetscViewer viewer;
1603ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1604ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1605ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1606ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1607ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1608ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1609ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1610ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1611ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1612ccaff030SJeremy L Thompson   }
161384d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1614ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1615ccaff030SJeremy L Thompson   }
1616ccaff030SJeremy L Thompson 
1617ccaff030SJeremy L Thompson   // Solve
1618ccaff030SJeremy L Thompson   start = MPI_Wtime();
1619ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1620ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1621ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1622ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1623ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1624ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
162584d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1626ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
162784d34d69SLeila Ghaffari                        "Time taken for solution (sec): %g\n",
1628ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1629ccaff030SJeremy L Thompson   }
1630ccaff030SJeremy L Thompson 
1631ccaff030SJeremy L Thompson   // Get error
163284d34d69SLeila Ghaffari   if (problem->non_zero_time && testChoice == TEST_NONE) {
1633ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1634ccaff030SJeremy L Thompson     PetscReal norm;
1635ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1636ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1637ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1638ccaff030SJeremy L Thompson 
163984d34d69SLeila Ghaffari     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
164084d34d69SLeila Ghaffari                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1641ccaff030SJeremy L Thompson 
1642ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1643ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1644cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1645ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1646ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1647ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
164884d34d69SLeila Ghaffari     // Clean up vectors
164984d34d69SLeila Ghaffari     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
165084d34d69SLeila Ghaffari     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1651ccaff030SJeremy L Thompson   }
1652ccaff030SJeremy L Thompson 
1653ccaff030SJeremy L Thompson   // Output Statistics
1654ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr);
165584d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1656ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1657ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1658ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1659ccaff030SJeremy L Thompson   }
166084d34d69SLeila Ghaffari   // Output numerical values from command line
166184d34d69SLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
166284d34d69SLeila Ghaffari 
166384d34d69SLeila Ghaffari   // compare reference solution values with current run
166484d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {
166584d34d69SLeila Ghaffari     PetscViewer viewer;
166684d34d69SLeila Ghaffari     // Read reference file
166784d34d69SLeila Ghaffari     Vec Qref;
166884d34d69SLeila Ghaffari     PetscReal error, Qrefnorm;
166984d34d69SLeila Ghaffari     ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
167084d34d69SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer);
167184d34d69SLeila Ghaffari     CHKERRQ(ierr);
167284d34d69SLeila Ghaffari     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
167384d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
167484d34d69SLeila Ghaffari 
167584d34d69SLeila Ghaffari     // Compute error with respect to reference solution
167684d34d69SLeila Ghaffari     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
167784d34d69SLeila Ghaffari     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
167884d34d69SLeila Ghaffari     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
167984d34d69SLeila Ghaffari     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
168084d34d69SLeila Ghaffari     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
168184d34d69SLeila Ghaffari     // Check error
168284d34d69SLeila Ghaffari     if (error > test->testtol) {
168384d34d69SLeila Ghaffari       ierr = PetscPrintf(PETSC_COMM_WORLD,
168484d34d69SLeila Ghaffari                          "Test failed with error norm %g\n",
168584d34d69SLeila Ghaffari                          (double)error); CHKERRQ(ierr);
168684d34d69SLeila Ghaffari     }
168784d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
168884d34d69SLeila Ghaffari   }
16899cf88b28Svaleriabarra 
1690ccaff030SJeremy L Thompson   // Clean up libCEED
16918b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1692ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1693ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1694ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1695ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
169684d34d69SLeila Ghaffari   CeedBasisDestroy(&basisq);
169784d34d69SLeila Ghaffari   CeedBasisDestroy(&basisx);
169884d34d69SLeila Ghaffari   CeedBasisDestroy(&basisxc);
169984d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictq);
170084d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictx);
170184d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdi);
1702ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1703ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1704ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1705ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1706ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1707ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
17086a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
17096a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
17106a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
17116a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
17126a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
17136a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
17146a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
1715*7659d40cSLeila Ghaffari   CeedQFunctionDestroy(&qf_Sur);
1716ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1717ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1718ccaff030SJeremy L Thompson 
1719ccaff030SJeremy L Thompson   // Clean up PETSc
1720ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1721ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1722ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1723ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1724ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1725ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1726ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1727ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1728ccaff030SJeremy L Thompson   return PetscFinalize();
1729ccaff030SJeremy L Thompson }
1730ccaff030SJeremy L Thompson 
1731