xref: /libCEED/examples/fluids/navierstokes.c (revision 9cf88b28e874a91fbe464b60f41b3973b01eeda4)
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 //
31ff6701fcSJed Brown //     ./navierstokes -ceed /cpu/self -problem density_current -degree 1
32ff6701fcSJed Brown //     ./navierstokes -ceed /gpu/occa -problem advection -degree 1
33ccaff030SJeremy L Thompson //
34*9cf88b28Svaleriabarra //TESTARGS -ceed {ceed_resource} -test -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
35*9cf88b28Svaleriabarra //TESTARGS -ceed {ceed_resource} -test -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
36*9cf88b28Svaleriabarra //TESTARGS -ceed {ceed_resource} -test -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"
49ccaff030SJeremy L Thompson #include "advection.h"
50ccaff030SJeremy L Thompson #include "advection2d.h"
51ccaff030SJeremy L Thompson #include "densitycurrent.h"
52ccaff030SJeremy L Thompson 
53ccaff030SJeremy L Thompson // Problem Options
54ccaff030SJeremy L Thompson typedef enum {
55ccaff030SJeremy L Thompson   NS_DENSITY_CURRENT = 0,
56ccaff030SJeremy L Thompson   NS_ADVECTION = 1,
57ccaff030SJeremy L Thompson   NS_ADVECTION2D = 2,
58ccaff030SJeremy L Thompson } problemType;
59ccaff030SJeremy L Thompson static const char *const problemTypes[] = {
60ccaff030SJeremy L Thompson   "density_current",
61ccaff030SJeremy L Thompson   "advection",
62ccaff030SJeremy L Thompson   "advection2d",
63ccaff030SJeremy L Thompson   "problemType","NS_",0
64ccaff030SJeremy L Thompson };
65ccaff030SJeremy L Thompson 
66ccaff030SJeremy L Thompson typedef enum {
67ccaff030SJeremy L Thompson   STAB_NONE = 0,
68ccaff030SJeremy L Thompson   STAB_SU = 1,   // Streamline Upwind
69ccaff030SJeremy L Thompson   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
70ccaff030SJeremy L Thompson } StabilizationType;
71ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = {
72ccaff030SJeremy L Thompson   "NONE",
73ccaff030SJeremy L Thompson   "SU",
74ccaff030SJeremy L Thompson   "SUPG",
75ccaff030SJeremy L Thompson   "StabilizationType", "STAB_", NULL
76ccaff030SJeremy L Thompson };
77ccaff030SJeremy L Thompson 
78ccaff030SJeremy L Thompson // Problem specific data
79ccaff030SJeremy L Thompson typedef struct {
80ccaff030SJeremy L Thompson   CeedInt dim, qdatasize;
81ccaff030SJeremy L Thompson   CeedQFunctionUser setup, ics, apply_rhs, apply_ifunction;
82ccaff030SJeremy L Thompson   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
83ccaff030SJeremy L Thompson                        PetscScalar[], void *);
84ccaff030SJeremy L Thompson   const char *setup_loc, *ics_loc, *apply_rhs_loc, *apply_ifunction_loc;
85ccaff030SJeremy L Thompson   const bool non_zero_time;
86ccaff030SJeremy L Thompson } problemData;
87ccaff030SJeremy L Thompson 
88ccaff030SJeremy L Thompson problemData problemOptions[] = {
89ccaff030SJeremy L Thompson   [NS_DENSITY_CURRENT] = {
90ccaff030SJeremy L Thompson     .dim                 = 3,
91ccaff030SJeremy L Thompson     .qdatasize           = 10,
92ccaff030SJeremy L Thompson     .setup               = Setup,
93ccaff030SJeremy L Thompson     .setup_loc           = Setup_loc,
94ccaff030SJeremy L Thompson     .ics                 = ICsDC,
95ccaff030SJeremy L Thompson     .ics_loc             = ICsDC_loc,
96ccaff030SJeremy L Thompson     .apply_rhs           = DC,
97ccaff030SJeremy L Thompson     .apply_rhs_loc       = DC_loc,
98ccaff030SJeremy L Thompson     .apply_ifunction     = IFunction_DC,
99ccaff030SJeremy L Thompson     .apply_ifunction_loc = IFunction_DC_loc,
100ccaff030SJeremy L Thompson     .bc                  = Exact_DC,
101ccaff030SJeremy L Thompson     .non_zero_time       = false,
102ccaff030SJeremy L Thompson   },
103ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
104ccaff030SJeremy L Thompson     .dim                 = 3,
105ccaff030SJeremy L Thompson     .qdatasize           = 10,
106ccaff030SJeremy L Thompson     .setup               = Setup,
107ccaff030SJeremy L Thompson     .setup_loc           = Setup_loc,
108ccaff030SJeremy L Thompson     .ics                 = ICsAdvection,
109ccaff030SJeremy L Thompson     .ics_loc             = ICsAdvection_loc,
110ccaff030SJeremy L Thompson     .apply_rhs           = Advection,
111ccaff030SJeremy L Thompson     .apply_rhs_loc       = Advection_loc,
112ccaff030SJeremy L Thompson     .apply_ifunction     = IFunction_Advection,
113ccaff030SJeremy L Thompson     .apply_ifunction_loc = IFunction_Advection_loc,
114ccaff030SJeremy L Thompson     .bc                  = Exact_Advection,
115ec9acdb8Svaleriabarra     .non_zero_time       = false,
116ccaff030SJeremy L Thompson   },
117ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
118ccaff030SJeremy L Thompson     .dim                 = 2,
119ccaff030SJeremy L Thompson     .qdatasize           = 5,
120ccaff030SJeremy L Thompson     .setup               = Setup2d,
121ccaff030SJeremy L Thompson     .setup_loc           = Setup2d_loc,
122ccaff030SJeremy L Thompson     .ics                 = ICsAdvection2d,
123ccaff030SJeremy L Thompson     .ics_loc             = ICsAdvection2d_loc,
124ccaff030SJeremy L Thompson     .apply_rhs           = Advection2d,
125ccaff030SJeremy L Thompson     .apply_rhs_loc       = Advection2d_loc,
126ccaff030SJeremy L Thompson     .apply_ifunction     = IFunction_Advection2d,
127ccaff030SJeremy L Thompson     .apply_ifunction_loc = IFunction_Advection2d_loc,
128ccaff030SJeremy L Thompson     .bc                  = Exact_Advection2d,
129ccaff030SJeremy L Thompson     .non_zero_time       = true,
130ccaff030SJeremy L Thompson   },
131ccaff030SJeremy L Thompson };
132ccaff030SJeremy L Thompson 
133ccaff030SJeremy L Thompson // PETSc user data
134ccaff030SJeremy L Thompson typedef struct User_ *User;
135ccaff030SJeremy L Thompson typedef struct Units_ *Units;
136ccaff030SJeremy L Thompson 
137ccaff030SJeremy L Thompson struct User_ {
138ccaff030SJeremy L Thompson   MPI_Comm comm;
139ccaff030SJeremy L Thompson   PetscInt outputfreq;
140ccaff030SJeremy L Thompson   DM dm;
141ccaff030SJeremy L Thompson   DM dmviz;
142ccaff030SJeremy L Thompson   Mat interpviz;
143ccaff030SJeremy L Thompson   Ceed ceed;
144ccaff030SJeremy L Thompson   Units units;
145ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
146ccaff030SJeremy L Thompson   CeedOperator op_rhs, op_ifunction;
147ccaff030SJeremy L Thompson   Vec M;
148ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
149ccaff030SJeremy L Thompson   PetscInt contsteps;
150ccaff030SJeremy L Thompson };
151ccaff030SJeremy L Thompson 
152ccaff030SJeremy L Thompson struct Units_ {
153ccaff030SJeremy L Thompson   // fundamental units
154ccaff030SJeremy L Thompson   PetscScalar meter;
155ccaff030SJeremy L Thompson   PetscScalar kilogram;
156ccaff030SJeremy L Thompson   PetscScalar second;
157ccaff030SJeremy L Thompson   PetscScalar Kelvin;
158ccaff030SJeremy L Thompson   // derived units
159ccaff030SJeremy L Thompson   PetscScalar Pascal;
160ccaff030SJeremy L Thompson   PetscScalar JperkgK;
161ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
162ccaff030SJeremy L Thompson   PetscScalar WpermK;
163ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
164ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
165ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
166ccaff030SJeremy L Thompson };
167ccaff030SJeremy L Thompson 
168ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
169ccaff030SJeremy L Thompson struct SimpleBC_ {
170ccaff030SJeremy L Thompson   PetscInt nwall, nslip[3];
171c063f476Svaleriabarra   PetscInt walls[6], slips[3][6];
17207af6069Svaleriabarra   PetscBool userbc;
173ccaff030SJeremy L Thompson };
174ccaff030SJeremy L Thompson 
175ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
176ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
177ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
178ccaff030SJeremy L Thompson }
179ccaff030SJeremy L Thompson 
180ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
181ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
182ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
183ccaff030SJeremy L Thompson 
184ccaff030SJeremy L Thompson   PetscSection   section;
185ccaff030SJeremy L Thompson   PetscInt       c, cStart, cEnd, Nelem, Ndof, *erestrict, eoffset, nfields, dim;
186ccaff030SJeremy L Thompson   PetscErrorCode ierr;
187ccaff030SJeremy L Thompson   Vec Uloc;
188ccaff030SJeremy L Thompson 
189ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
190ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
191ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm,&section); CHKERRQ(ierr);
192ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
193ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
194ccaff030SJeremy L Thompson   fieldoff[0] = 0;
195ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
196ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
197ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
198ccaff030SJeremy L Thompson   }
199ccaff030SJeremy L Thompson 
200ccaff030SJeremy L Thompson   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
201ccaff030SJeremy L Thompson   Nelem = cEnd - cStart;
202ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
203ccaff030SJeremy L Thompson   for (c=cStart,eoffset=0; c<cEnd; c++) {
204ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
205ccaff030SJeremy L Thompson     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
206ccaff030SJeremy L Thompson                                    &indices, NULL); CHKERRQ(ierr);
207ccaff030SJeremy L Thompson     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
208ccaff030SJeremy L Thompson           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
209ccaff030SJeremy L Thompson           c);
210ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
211ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
212ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
213ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
214ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
215ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
216ccaff030SJeremy L Thompson           if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j])
217ccaff030SJeremy L Thompson               != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j)
218ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
219ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
220ccaff030SJeremy L Thompson                      c, i, f, j);
221ccaff030SJeremy L Thompson         }
222ccaff030SJeremy L Thompson       }
223ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
224ccaff030SJeremy L Thompson       PetscInt loc = Involute(indices[i*ncomp[0]]);
225ccaff030SJeremy L Thompson       erestrict[eoffset++] = loc / fieldoff[nfields];
226ccaff030SJeremy L Thompson     }
227ccaff030SJeremy L Thompson     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
228ccaff030SJeremy L Thompson                                        &indices, NULL); CHKERRQ(ierr);
229ccaff030SJeremy L Thompson   }
230ccaff030SJeremy L Thompson   if (eoffset != Nelem*PetscPowInt(P, dim)) SETERRQ3(PETSC_COMM_SELF,
231ccaff030SJeremy L Thompson         PETSC_ERR_LIB, "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
232ccaff030SJeremy L Thompson         PetscPowInt(P, dim),eoffset);
233ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
234ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
235ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
236ccaff030SJeremy L Thompson   CeedElemRestrictionCreate(ceed, CEED_INTERLACED, Nelem, PetscPowInt(P, dim),
237ccaff030SJeremy L Thompson                             Ndof/fieldoff[nfields], fieldoff[nfields],
238ccaff030SJeremy L Thompson                             CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, Erestrict);
239ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
240ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
241ccaff030SJeremy L Thompson }
242ccaff030SJeremy L Thompson 
243ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
244ccaff030SJeremy L Thompson   PetscErrorCode ierr;
245ccaff030SJeremy L Thompson   PetscInt m;
246ccaff030SJeremy L Thompson 
247ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
248ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
249ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
250ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
251ccaff030SJeremy L Thompson }
252ccaff030SJeremy L Thompson 
253ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
254ccaff030SJeremy L Thompson   PetscErrorCode ierr;
255ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
256ccaff030SJeremy L Thompson   PetscScalar *a;
257ccaff030SJeremy L Thompson 
258ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
259ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
260ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
261ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
262380e0a58Svaleriabarra                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
263380e0a58Svaleriabarra                                   mpetsc, mceed);
264ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
265ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
266ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
267ccaff030SJeremy L Thompson }
268ccaff030SJeremy L Thompson 
269ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
270ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
271ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
272ccaff030SJeremy L Thompson   PetscErrorCode ierr;
273ccaff030SJeremy L Thompson   Vec Qbc;
274ccaff030SJeremy L Thompson 
275ccaff030SJeremy L Thompson   PetscFunctionBegin;
276ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
277ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
278ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
279ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
280ccaff030SJeremy L Thompson }
281ccaff030SJeremy L Thompson 
282ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
283ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
284ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
285ccaff030SJeremy L Thompson   PetscErrorCode ierr;
286ccaff030SJeremy L Thompson   User user = *(User *)userData;
287ccaff030SJeremy L Thompson   PetscScalar *q, *g;
288ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
289ccaff030SJeremy L Thompson 
290ccaff030SJeremy L Thompson   // Global-to-local
291ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
292ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
293ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
294ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
295ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
296ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
297ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
298ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
299ccaff030SJeremy L Thompson 
300ccaff030SJeremy L Thompson   // Ceed Vectors
301ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
302ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
303ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
304ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
305ccaff030SJeremy L Thompson 
306ccaff030SJeremy L Thompson   // Apply CEED operator
307ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
308ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
309ccaff030SJeremy L Thompson 
310ccaff030SJeremy L Thompson   // Restore vectors
311ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
312ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
313ccaff030SJeremy L Thompson 
314ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
315ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
316ccaff030SJeremy L Thompson 
317ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
318ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
319ccaff030SJeremy L Thompson   CHKERRQ(ierr);
320ccaff030SJeremy L Thompson 
321ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
322ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
323ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
324ccaff030SJeremy L Thompson }
325ccaff030SJeremy L Thompson 
326ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
327ccaff030SJeremy L Thompson                                    void *userData) {
328ccaff030SJeremy L Thompson   PetscErrorCode ierr;
329ccaff030SJeremy L Thompson   User user = *(User *)userData;
330ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
331ccaff030SJeremy L Thompson   PetscScalar *g;
332ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
333ccaff030SJeremy L Thompson 
334ccaff030SJeremy L Thompson   // Global-to-local
335ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
336ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
337ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
338ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
339ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
340ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
341ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
342ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
343ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
344ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
345ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
346ccaff030SJeremy L Thompson 
347ccaff030SJeremy L Thompson   // Ceed Vectors
348ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
349ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
350ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
351ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
352ccaff030SJeremy L Thompson                      (PetscScalar *)q);
353ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
354ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
355ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
356ccaff030SJeremy L Thompson 
357ccaff030SJeremy L Thompson   // Apply CEED operator
358ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
359ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
360ccaff030SJeremy L Thompson 
361ccaff030SJeremy L Thompson   // Restore vectors
362ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
363ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
364ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
365ccaff030SJeremy L Thompson 
366ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
367ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
368ccaff030SJeremy L Thompson 
369ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
370ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
371ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
372ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
373ccaff030SJeremy L Thompson }
374ccaff030SJeremy L Thompson 
375ccaff030SJeremy L Thompson // User provided TS Monitor
376ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
377ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
378ccaff030SJeremy L Thompson   User user = ctx;
379ccaff030SJeremy L Thompson   Vec Qloc;
380ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
381ccaff030SJeremy L Thompson   PetscViewer viewer;
382ccaff030SJeremy L Thompson   PetscErrorCode ierr;
383ccaff030SJeremy L Thompson 
384ccaff030SJeremy L Thompson   // Set up output
385ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
386ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
387ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
388ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
389ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
390ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
391ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
392ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
393ccaff030SJeremy L Thompson 
394ccaff030SJeremy L Thompson   // Output
395ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
396ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
397ccaff030SJeremy L Thompson   CHKERRQ(ierr);
398ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
399ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
400ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
401ccaff030SJeremy L Thompson   if (user->dmviz) {
402ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
403ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
404ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
405ccaff030SJeremy L Thompson 
406ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
407ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
408ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
409ccaff030SJeremy L Thompson     CHKERRQ(ierr);
410ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
411ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
412ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
413ccaff030SJeremy L Thompson     CHKERRQ(ierr);
414ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
415ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
416ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
417ccaff030SJeremy L Thompson     CHKERRQ(ierr);
418ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
419ccaff030SJeremy L Thompson                               filepath_refined,
420ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
421ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
422ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
423ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
424ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
425ccaff030SJeremy L Thompson   }
426ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
427ccaff030SJeremy L Thompson 
428ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
429ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
430ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
431ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
432ccaff030SJeremy L Thompson   CHKERRQ(ierr);
433ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
434ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
435ccaff030SJeremy L Thompson 
436ccaff030SJeremy L Thompson   // Save time stamp
437ccaff030SJeremy L Thompson   // Dimensionalize time back
438ccaff030SJeremy L Thompson   time /= user->units->second;
439ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
440ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
441ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
442ccaff030SJeremy L Thompson   CHKERRQ(ierr);
443ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
444ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
445ccaff030SJeremy L Thompson   #else
446ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
447ccaff030SJeremy L Thompson   #endif
448ccaff030SJeremy L Thompson   CHKERRQ(ierr);
449ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
450ccaff030SJeremy L Thompson 
451ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
452ccaff030SJeremy L Thompson }
453ccaff030SJeremy L Thompson 
45427dd567eSvaleriabarra static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
455ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
456ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
457ccaff030SJeremy L Thompson   PetscErrorCode ierr;
458ccaff030SJeremy L Thompson   CeedVector multlvec;
459ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
460ccaff030SJeremy L Thompson 
461ccaff030SJeremy L Thompson   ctxSetup->time = time;
462ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
463ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
464ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
465ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
466ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
467ccaff030SJeremy L Thompson 
468ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
469ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
470ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
471ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
472ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
473ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
474ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
475ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
476ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
477ccaff030SJeremy L Thompson   CHKERRQ(ierr);
478ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
479ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
480ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
481ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
482ccaff030SJeremy L Thompson 
483ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
484ccaff030SJeremy L Thompson }
485ccaff030SJeremy L Thompson 
486ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
487ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
488ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
489ccaff030SJeremy L Thompson   PetscErrorCode ierr;
490ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
491ccaff030SJeremy L Thompson   CeedOperator op_mass;
492ccaff030SJeremy L Thompson   CeedVector mceed;
493ccaff030SJeremy L Thompson   Vec Mloc;
494ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
495ccaff030SJeremy L Thompson 
496ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
497ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
498ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
499ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
500ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
501ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
502ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
503ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
504ccaff030SJeremy L Thompson 
505ccaff030SJeremy L Thompson   // Create the mass operator
506ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
507ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
508ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
509ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
510ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
511ccaff030SJeremy L Thompson 
512ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
513ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
514ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
515ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
516ccaff030SJeremy L Thompson 
517ccaff030SJeremy L Thompson   {
518ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
519ccaff030SJeremy L Thompson     CeedVector onesvec;
520ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
521ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
522ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
523ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
524ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
525ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
526ccaff030SJeremy L Thompson   }
527ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
528ccaff030SJeremy L Thompson 
529ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
530ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
531ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
532ccaff030SJeremy L Thompson 
533ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
534ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
535ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
536ccaff030SJeremy L Thompson }
537ccaff030SJeremy L Thompson 
538c063f476Svaleriabarra static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
539ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
540ccaff030SJeremy L Thompson   PetscErrorCode ierr;
541ccaff030SJeremy L Thompson 
542ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
543ccaff030SJeremy L Thompson   {
544ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
545ccaff030SJeremy L Thompson     PetscFE fe;
546ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
547ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
548ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
549ff6701fcSJed Brown                                  &fe);
550ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
551ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
552ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
55307af6069Svaleriabarra     {
55407af6069Svaleriabarra       PetscInt comps[1] = {1};
55507af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
55607af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
55707af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
55807af6069Svaleriabarra       comps[0] = 2;
55907af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
56007af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
56107af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
56207af6069Svaleriabarra       comps[0] = 3;
56307af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
56407af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
56507af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
56607af6069Svaleriabarra     }
56707af6069Svaleriabarra     if (bc->userbc == PETSC_TRUE) {
56807af6069Svaleriabarra       for (PetscInt c = 0; c < 3; c++) {
56907af6069Svaleriabarra         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
57007af6069Svaleriabarra           for (PetscInt w = 0; w < bc->nwall; w++) {
57107af6069Svaleriabarra             if (bc->slips[c][s] == bc->walls[w])
57207af6069Svaleriabarra               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
57307af6069Svaleriabarra                        "Boundary condition already set on face %D!\n", bc->walls[w]);
57407af6069Svaleriabarra 
57507af6069Svaleriabarra           }
57607af6069Svaleriabarra         }
57707af6069Svaleriabarra       }
57807af6069Svaleriabarra     }
57935d2ed1bSvaleriabarra     // Wall boundary conditions are zero energy density and zero flux for
58035d2ed1bSvaleriabarra     //   velocity in advection/advection2d, and zero velocity and zero flux
58135d2ed1bSvaleriabarra     //   for mass density and energy density in density_current
582ccaff030SJeremy L Thompson     {
58335d2ed1bSvaleriabarra       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
58435d2ed1bSvaleriabarra         PetscInt comps[1] = {4};
58535d2ed1bSvaleriabarra         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
58635d2ed1bSvaleriabarra                              1, comps, (void(*)(void))problem->bc,
58735d2ed1bSvaleriabarra                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
58835d2ed1bSvaleriabarra       } else if (problem->bc == Exact_DC) {
589ccaff030SJeremy L Thompson         PetscInt comps[3] = {1, 2, 3};
590ccaff030SJeremy L Thompson         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
591ccaff030SJeremy L Thompson                              3, comps, (void(*)(void))problem->bc,
592ccaff030SJeremy L Thompson                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
59307af6069Svaleriabarra       } else
59407af6069Svaleriabarra         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
595c063f476Svaleriabarra                 "Undefined boundary conditions for this problem");
596ccaff030SJeremy L Thompson     }
597ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
598ccaff030SJeremy L Thompson     CHKERRQ(ierr);
599ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
600ccaff030SJeremy L Thompson   }
601ccaff030SJeremy L Thompson   {
602ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
603ccaff030SJeremy L Thompson     PetscSection section;
604ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
605ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
606ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
607ccaff030SJeremy L Thompson     CHKERRQ(ierr);
608ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
609ccaff030SJeremy L Thompson     CHKERRQ(ierr);
610ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
611ccaff030SJeremy L Thompson     CHKERRQ(ierr);
612ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
613ccaff030SJeremy L Thompson     CHKERRQ(ierr);
614ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
615ccaff030SJeremy L Thompson     CHKERRQ(ierr);
616ccaff030SJeremy L Thompson   }
617ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
618ccaff030SJeremy L Thompson }
619ccaff030SJeremy L Thompson 
620ccaff030SJeremy L Thompson int main(int argc, char **argv) {
621ccaff030SJeremy L Thompson   PetscInt ierr;
622ccaff030SJeremy L Thompson   MPI_Comm comm;
623380e0a58Svaleriabarra   DM dm, dmcoord, dmviz;
624ccaff030SJeremy L Thompson   Mat interpviz;
625ccaff030SJeremy L Thompson   TS ts;
626ccaff030SJeremy L Thompson   TSAdapt adapt;
627ccaff030SJeremy L Thompson   User user;
628ccaff030SJeremy L Thompson   Units units;
629ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
630ccaff030SJeremy L Thompson   PetscInt cStart, cEnd, localNelem, lnodes, steps;
631ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
632ccaff030SJeremy L Thompson   PetscMPIInt rank;
633ccaff030SJeremy L Thompson   PetscScalar ftime;
634ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
635ccaff030SJeremy L Thompson   Ceed ceed;
636ccaff030SJeremy L Thompson   CeedInt numP, numQ;
637ccaff030SJeremy L Thompson   CeedVector xcorners, qdata, q0ceed;
638ccaff030SJeremy L Thompson   CeedBasis basisx, basisxc, basisq;
639ccaff030SJeremy L Thompson   CeedElemRestriction restrictx, restrictxcoord, restrictq, restrictqdi;
640ccaff030SJeremy L Thompson   CeedQFunction qf_setup, qf_ics, qf_rhs, qf_ifunction;
641ccaff030SJeremy L Thompson   CeedOperator op_setup, op_ics;
642ccaff030SJeremy L Thompson   CeedScalar Rd;
643ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
644ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
645ccaff030SJeremy L Thompson   problemType problemChoice;
646ccaff030SJeremy L Thompson   problemData *problem = NULL;
647ccaff030SJeremy L Thompson   StabilizationType stab;
648cb3e2689Svaleriabarra   PetscBool   test, implicit;
649cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
650ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
65107af6069Svaleriabarra     .nslip = {2, 2, 2},
65207af6069Svaleriabarra     .slips = {{5, 6}, {3, 4}, {1, 2}}
653ccaff030SJeremy L Thompson   };
654ccaff030SJeremy L Thompson   double start, cpu_time_used;
655ccaff030SJeremy L Thompson 
656ccaff030SJeremy L Thompson   // Create the libCEED contexts
657ccaff030SJeremy L Thompson   PetscScalar meter     = 1e-2;     // 1 meter in scaled length units
658ccaff030SJeremy L Thompson   PetscScalar second    = 1e-2;     // 1 second in scaled time units
659ccaff030SJeremy L Thompson   PetscScalar kilogram  = 1e-6;     // 1 kilogram in scaled mass units
660ccaff030SJeremy L Thompson   PetscScalar Kelvin    = 1;        // 1 Kelvin in scaled temperature units
661ccaff030SJeremy L Thompson   CeedScalar theta0     = 300.;     // K
662ccaff030SJeremy L Thompson   CeedScalar thetaC     = -15.;     // K
663ccaff030SJeremy L Thompson   CeedScalar P0         = 1.e5;     // Pa
664ccaff030SJeremy L Thompson   CeedScalar N          = 0.01;     // 1/s
665ccaff030SJeremy L Thompson   CeedScalar cv         = 717.;     // J/(kg K)
666ccaff030SJeremy L Thompson   CeedScalar cp         = 1004.;    // J/(kg K)
667ccaff030SJeremy L Thompson   CeedScalar g          = 9.81;     // m/s^2
668ccaff030SJeremy L Thompson   CeedScalar lambda     = -2./3.;   // -
669ccaff030SJeremy L Thompson   CeedScalar mu         = 75.;      // Pa s, dynamic viscosity
670ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
671ccaff030SJeremy L Thompson   CeedScalar k          = 0.02638;  // W/(m K)
672ccaff030SJeremy L Thompson   CeedScalar CtauS      = 0.;       // dimensionless
673ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
674ccaff030SJeremy L Thompson   PetscScalar lx        = 8000.;    // m
675ccaff030SJeremy L Thompson   PetscScalar ly        = 8000.;    // m
676ccaff030SJeremy L Thompson   PetscScalar lz        = 4000.;    // m
677ccaff030SJeremy L Thompson   CeedScalar rc         = 1000.;    // m (Radius of bubble)
678ccaff030SJeremy L Thompson   PetscScalar resx      = 1000.;    // m (resolution in x)
679ccaff030SJeremy L Thompson   PetscScalar resy      = 1000.;    // m (resolution in y)
680ccaff030SJeremy L Thompson   PetscScalar resz      = 1000.;    // m (resolution in z)
681ccaff030SJeremy L Thompson   PetscInt outputfreq   = 10;       // -
682ccaff030SJeremy L Thompson   PetscInt contsteps    = 0;        // -
683ff6701fcSJed Brown   PetscInt degree       = 1;        // -
684ccaff030SJeremy L Thompson   PetscInt qextra       = 2;        // -
685ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
686ccaff030SJeremy L Thompson 
687ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
688ccaff030SJeremy L Thompson   if (ierr) return ierr;
689ccaff030SJeremy L Thompson 
690ccaff030SJeremy L Thompson   // Allocate PETSc context
691ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
692ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
693ccaff030SJeremy L Thompson 
694ccaff030SJeremy L Thompson   // Parse command line options
695ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
696ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
697ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
698ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
699ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
700ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
701ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-test", "Run in test mode",
702ccaff030SJeremy L Thompson                           NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr);
703ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
704ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
705ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
706ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
707ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
708ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
709ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
710ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
711ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
712ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
713ccaff030SJeremy L Thompson   CHKERRQ(ierr);
7148c662231Svaleriabarra   if (!implicit && stab != STAB_NONE) {
7158c662231Svaleriabarra     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
7168c662231Svaleriabarra     CHKERRQ(ierr);
7178c662231Svaleriabarra   }
718ccaff030SJeremy L Thompson   {
719ccaff030SJeremy L Thompson     PetscInt len;
720ccaff030SJeremy L Thompson     PetscBool flg;
721ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
722ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
723ccaff030SJeremy L Thompson                                 NULL, bc.walls,
724ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
725ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
726ccaff030SJeremy L Thompson     if (flg) bc.nwall = len;
727ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
728ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
729ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
730ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
731ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
732ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
733ccaff030SJeremy L Thompson                                    &len), &flg);
734ccaff030SJeremy L Thompson       CHKERRQ(ierr);
73507af6069Svaleriabarra       if (flg) {
73607af6069Svaleriabarra         bc.nslip[j] = len;
73707af6069Svaleriabarra         bc.userbc = PETSC_TRUE;
73807af6069Svaleriabarra       }
739ccaff030SJeremy L Thompson     }
740ccaff030SJeremy L Thompson   }
741cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
742cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
743cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
744ccaff030SJeremy L Thompson   CHKERRQ(ierr);
745ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
746ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
747ccaff030SJeremy L Thompson   meter = fabs(meter);
748ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
749ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
750ccaff030SJeremy L Thompson   second = fabs(second);
751ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
752ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
753ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
754ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
755ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
756ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
757ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
758ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
759ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
760ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
761ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
762ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
763ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
764ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
765ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
766ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
767ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
768ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
769ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
770ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
771ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
772ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
773ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
774ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
775ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
776ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
777ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
778ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
779ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
780ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
781ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
78240c555b9Svaleriabarra   if (stab == STAB_NONE && CtauS != 0) {
783c063f476Svaleriabarra     ierr = PetscPrintf(comm,
784c063f476Svaleriabarra                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
78540c555b9Svaleriabarra     CHKERRQ(ierr);
78640c555b9Svaleriabarra   }
787ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
788ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
789ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
790ccaff030SJeremy L Thompson   CHKERRQ(ierr);
791ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
792ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
793ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
794ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
795ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
796ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
797ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
798ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
799ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
800ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
801ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
802ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
803ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
804ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
805ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
806ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
807ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
808ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
809ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
810ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
811ccaff030SJeremy L Thompson   n = problem->dim;
812ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
813ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
814ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
815ccaff030SJeremy L Thompson   {
816ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
817ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
818ccaff030SJeremy L Thompson     if (norm > 0) {
819ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
820ccaff030SJeremy L Thompson     }
821ccaff030SJeremy L Thompson   }
822ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
823ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
824ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
825ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
826ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
827ff6701fcSJed Brown   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
828ff6701fcSJed Brown                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
829ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
830ccaff030SJeremy L Thompson                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
831ccaff030SJeremy L Thompson   PetscStrncpy(user->outputfolder, ".", 2);
832ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
833ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
834ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
835ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
836ccaff030SJeremy L Thompson 
837ccaff030SJeremy L Thompson   // Define derived units
838ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
839ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
840ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
841ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
842ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
843ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
844ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
845ccaff030SJeremy L Thompson 
846ccaff030SJeremy L Thompson   // Scale variables to desired units
847ccaff030SJeremy L Thompson   theta0 *= Kelvin;
848ccaff030SJeremy L Thompson   thetaC *= Kelvin;
849ccaff030SJeremy L Thompson   P0 *= Pascal;
850ccaff030SJeremy L Thompson   N *= (1./second);
851ccaff030SJeremy L Thompson   cv *= JperkgK;
852ccaff030SJeremy L Thompson   cp *= JperkgK;
853ccaff030SJeremy L Thompson   Rd = cp - cv;
854ccaff030SJeremy L Thompson   g *= mpersquareds;
855ccaff030SJeremy L Thompson   mu *= Pascal * second;
856ccaff030SJeremy L Thompson   k *= WpermK;
857ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
858ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
859ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
860ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
861ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
862ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
863ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
864ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
865ccaff030SJeremy L Thompson 
866ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
867ccaff030SJeremy L Thompson                 qdatasize = problem->qdatasize;
868ccaff030SJeremy L Thompson   // Set up the libCEED context
869ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
870ccaff030SJeremy L Thompson     .theta0 = theta0,
871ccaff030SJeremy L Thompson     .thetaC = thetaC,
872ccaff030SJeremy L Thompson     .P0 = P0,
873ccaff030SJeremy L Thompson     .N = N,
874ccaff030SJeremy L Thompson     .cv = cv,
875ccaff030SJeremy L Thompson     .cp = cp,
876ccaff030SJeremy L Thompson     .Rd = Rd,
877ccaff030SJeremy L Thompson     .g = g,
878ccaff030SJeremy L Thompson     .rc = rc,
879ccaff030SJeremy L Thompson     .lx = lx,
880ccaff030SJeremy L Thompson     .ly = ly,
881ccaff030SJeremy L Thompson     .lz = lz,
882ccaff030SJeremy L Thompson     .center[0] = center[0],
883ccaff030SJeremy L Thompson     .center[1] = center[1],
884ccaff030SJeremy L Thompson     .center[2] = center[2],
885ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
886ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
887ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
888ccaff030SJeremy L Thompson     .time = 0,
889ccaff030SJeremy L Thompson   };
890ccaff030SJeremy L Thompson 
891380e0a58Svaleriabarra   // Create the mesh
892ccaff030SJeremy L Thompson   {
893ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
894ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
895ff62f740Svaleriabarra                                NULL, PETSC_TRUE, &dm);
896ccaff030SJeremy L Thompson     CHKERRQ(ierr);
897ccaff030SJeremy L Thompson   }
898380e0a58Svaleriabarra 
899380e0a58Svaleriabarra   // Distribute the mesh over processes
900380e0a58Svaleriabarra   {
901ccaff030SJeremy L Thompson     DM               dmDist = NULL;
902ccaff030SJeremy L Thompson     PetscPartitioner part;
903ccaff030SJeremy L Thompson 
904ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
905ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
906ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
907ccaff030SJeremy L Thompson     if (dmDist) {
908ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
909ccaff030SJeremy L Thompson       dm  = dmDist;
910ccaff030SJeremy L Thompson     }
911ccaff030SJeremy L Thompson   }
912ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
913ccaff030SJeremy L Thompson 
914380e0a58Svaleriabarra   // Setup DM
915ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
916ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
917ff6701fcSJed Brown   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
918380e0a58Svaleriabarra 
919380e0a58Svaleriabarra   // Print FEM space information
920ccaff030SJeremy L Thompson   if (!test) {
921ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
922380e0a58Svaleriabarra                        "Degree of FEM space: %D\n",
923ff6701fcSJed Brown                        degree); CHKERRQ(ierr);
924ccaff030SJeremy L Thompson   }
925380e0a58Svaleriabarra 
926380e0a58Svaleriabarra   // Refine DM for high-order viz
927ccaff030SJeremy L Thompson   dmviz = NULL;
928ccaff030SJeremy L Thompson   interpviz = NULL;
929ccaff030SJeremy L Thompson   if (viz_refine) {
930ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
931ff6701fcSJed Brown 
932ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
933ff6701fcSJed Brown     dmhierarchy[0] = dm;
934ff6701fcSJed Brown     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
935ff6701fcSJed Brown       Mat interp_next;
936ff6701fcSJed Brown 
937ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
938ccaff030SJeremy L Thompson       CHKERRQ(ierr);
939ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
940ff6701fcSJed Brown       d = (d + 1) / 2;
941ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
942ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
943ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
944ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
945ff6701fcSJed Brown       if (!i) interpviz = interp_next;
946ff6701fcSJed Brown       else {
947ff6701fcSJed Brown         Mat C;
948ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
949ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
950ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
951ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
952ff6701fcSJed Brown         interpviz = C;
953ff6701fcSJed Brown       }
954ff6701fcSJed Brown     }
955cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
956ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
957cb3e2689Svaleriabarra     }
958ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
959ccaff030SJeremy L Thompson   }
960ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
961ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
962ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
963ccaff030SJeremy L Thompson   lnodes /= ncompq;
964ccaff030SJeremy L Thompson 
965ccaff030SJeremy L Thompson   {
966ccaff030SJeremy L Thompson     // Print grid information
967ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
968ccaff030SJeremy L Thompson     int comm_size;
969ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
970ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
971ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
972ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
973ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
974ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
975ccaff030SJeremy L Thompson     if (!test) {
976380e0a58Svaleriabarra       ierr = PetscPrintf(comm, "Global FEM dofs: %D (%D owned) on %d rank(s)\n",
977ccaff030SJeremy L Thompson                          gdofs, odofs, comm_size); CHKERRQ(ierr);
978ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm, "Local FEM nodes: %D\n", lnodes); CHKERRQ(ierr);
979ccaff030SJeremy L Thompson       ierr = PetscPrintf(comm, "dm_plex_box_faces: %s\n", box_faces_str);
980ccaff030SJeremy L Thompson       CHKERRQ(ierr);
981ccaff030SJeremy L Thompson     }
982ccaff030SJeremy L Thompson 
983ccaff030SJeremy L Thompson   }
984ccaff030SJeremy L Thompson 
985ccaff030SJeremy L Thompson   // Set up global mass vector
986ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
987ccaff030SJeremy L Thompson 
988ccaff030SJeremy L Thompson   // Set up CEED
989ccaff030SJeremy L Thompson   // CEED Bases
990ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
991ccaff030SJeremy L Thompson   numP = degree + 1;
992ccaff030SJeremy L Thompson   numQ = numP + qextra;
993ccaff030SJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
994ccaff030SJeremy L Thompson                                   &basisq);
995ccaff030SJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
996ccaff030SJeremy L Thompson                                   &basisx);
997ccaff030SJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
998ccaff030SJeremy L Thompson                                   CEED_GAUSS_LOBATTO, &basisxc);
999ccaff030SJeremy L Thompson 
1000ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1001ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1002ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1003ccaff030SJeremy L Thompson 
1004ccaff030SJeremy L Thompson   // CEED Restrictions
1005ccaff030SJeremy L Thompson   ierr = CreateRestrictionFromPlex(ceed, dm, degree+1, &restrictq);
1006ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1007ccaff030SJeremy L Thompson   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, &restrictx); CHKERRQ(ierr);
1008ccaff030SJeremy L Thompson   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
1009ccaff030SJeremy L Thompson   localNelem = cEnd - cStart;
1010ccaff030SJeremy L Thompson   CeedInt numQdim = CeedIntPow(numQ, dim);
1011ccaff030SJeremy L Thompson   CeedElemRestrictionCreateStrided(ceed, localNelem, numQdim,
1012ccaff030SJeremy L Thompson                                    localNelem*numQdim, qdatasize,
1013ccaff030SJeremy L Thompson                                    CEED_STRIDES_BACKEND, &restrictqdi);
1014ccaff030SJeremy L Thompson   CeedElemRestrictionCreateStrided(ceed, localNelem, PetscPowInt(numP, dim),
1015ccaff030SJeremy L Thompson                                    localNelem*PetscPowInt(numP, dim), ncompx,
1016ccaff030SJeremy L Thompson                                    CEED_STRIDES_BACKEND, &restrictxcoord);
1017ccaff030SJeremy L Thompson 
1018ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1019ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1020ccaff030SJeremy L Thompson 
1021ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1022ccaff030SJeremy L Thompson   CeedInt Nqpts;
1023ccaff030SJeremy L Thompson   CeedBasisGetNumQuadraturePoints(basisq, &Nqpts);
1024ccaff030SJeremy L Thompson   CeedVectorCreate(ceed, qdatasize*localNelem*Nqpts, &qdata);
1025ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1026ccaff030SJeremy L Thompson 
1027ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1028ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->setup, problem->setup_loc,
1029ccaff030SJeremy L Thompson                               &qf_setup);
1030ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "dx", ncompx*dim, CEED_EVAL_GRAD);
1031ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT);
1032ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup, "qdata", qdatasize, CEED_EVAL_NONE);
1033ccaff030SJeremy L Thompson 
1034ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1035ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1036ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1037ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1038ccaff030SJeremy L Thompson 
1039ccaff030SJeremy L Thompson   qf_rhs = NULL;
1040ccaff030SJeremy L Thompson   if (problem->apply_rhs) { // Create the Q-function that defines the action of the RHS operator
1041ccaff030SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_rhs,
1042ccaff030SJeremy L Thompson                                 problem->apply_rhs_loc, &qf_rhs);
1043ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_rhs, "q", ncompq, CEED_EVAL_INTERP);
1044ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_rhs, "dq", ncompq*dim, CEED_EVAL_GRAD);
1045ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_rhs, "qdata", qdatasize, CEED_EVAL_NONE);
1046ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_rhs, "x", ncompx, CEED_EVAL_INTERP);
1047ccaff030SJeremy L Thompson     CeedQFunctionAddOutput(qf_rhs, "v", ncompq, CEED_EVAL_INTERP);
1048ccaff030SJeremy L Thompson     CeedQFunctionAddOutput(qf_rhs, "dv", ncompq*dim, CEED_EVAL_GRAD);
1049ccaff030SJeremy L Thompson   }
1050ccaff030SJeremy L Thompson 
1051ccaff030SJeremy L Thompson   qf_ifunction = NULL;
1052ccaff030SJeremy L Thompson   if (problem->apply_ifunction) { // Create the Q-function that defines the action of the IFunction
1053ccaff030SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_ifunction,
1054ccaff030SJeremy L Thompson                                 problem->apply_ifunction_loc, &qf_ifunction);
1055ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_ifunction, "q", ncompq, CEED_EVAL_INTERP);
1056ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_ifunction, "dq", ncompq*dim, CEED_EVAL_GRAD);
1057ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_ifunction, "qdot", ncompq, CEED_EVAL_INTERP);
1058ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_ifunction, "qdata", qdatasize, CEED_EVAL_NONE);
1059ccaff030SJeremy L Thompson     CeedQFunctionAddInput(qf_ifunction, "x", ncompx, CEED_EVAL_INTERP);
1060ccaff030SJeremy L Thompson     CeedQFunctionAddOutput(qf_ifunction, "v", ncompq, CEED_EVAL_INTERP);
1061ccaff030SJeremy L Thompson     CeedQFunctionAddOutput(qf_ifunction, "dv", ncompq*dim, CEED_EVAL_GRAD);
1062ccaff030SJeremy L Thompson   }
1063ccaff030SJeremy L Thompson 
1064ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1065ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup);
1066ccaff030SJeremy L Thompson   CeedOperatorSetField(op_setup, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1067ccaff030SJeremy L Thompson   CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE,
1068ccaff030SJeremy L Thompson                        basisx, CEED_VECTOR_NONE);
1069ccaff030SJeremy L Thompson   CeedOperatorSetField(op_setup, "qdata", restrictqdi,
1070ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1071ccaff030SJeremy L Thompson 
1072ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1073ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
1074ccaff030SJeremy L Thompson   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
1075ccaff030SJeremy L Thompson   CeedOperatorSetField(op_ics, "q0", restrictq,
1076ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1077ccaff030SJeremy L Thompson 
1078ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
1079ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
1080ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1081ccaff030SJeremy L Thompson 
1082ccaff030SJeremy L Thompson   if (qf_rhs) { // Create the RHS physics operator
1083ccaff030SJeremy L Thompson     CeedOperator op;
1084ccaff030SJeremy L Thompson     CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op);
1085ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
1086ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
1087ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "qdata", restrictqdi,
1088ccaff030SJeremy L Thompson                          CEED_BASIS_COLLOCATED, qdata);
1089ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
1090ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
1091ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1092ccaff030SJeremy L Thompson     user->op_rhs = op;
1093ccaff030SJeremy L Thompson   }
1094ccaff030SJeremy L Thompson 
1095ccaff030SJeremy L Thompson   if (qf_ifunction) { // Create the IFunction operator
1096ccaff030SJeremy L Thompson     CeedOperator op;
1097ccaff030SJeremy L Thompson     CeedOperatorCreate(ceed, qf_ifunction, NULL, NULL, &op);
1098ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
1099ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
1100ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
1101ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "qdata", restrictqdi,
1102ccaff030SJeremy L Thompson                          CEED_BASIS_COLLOCATED, qdata);
1103ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
1104ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
1105ccaff030SJeremy L Thompson     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1106ccaff030SJeremy L Thompson     user->op_ifunction = op;
1107ccaff030SJeremy L Thompson   }
1108ccaff030SJeremy L Thompson 
1109ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1110ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1111ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1112ccaff030SJeremy L Thompson     .CtauS = CtauS,
1113ccaff030SJeremy L Thompson     .strong_form = strong_form,
1114ccaff030SJeremy L Thompson     .stabilization = stab,
1115ccaff030SJeremy L Thompson   };
1116ccaff030SJeremy L Thompson   switch (problemChoice) {
1117ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1118ccaff030SJeremy L Thompson     if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxNS, sizeof ctxNS);
1119ccaff030SJeremy L Thompson     if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxNS,
1120ccaff030SJeremy L Thompson           sizeof ctxNS);
1121ccaff030SJeremy L Thompson     break;
1122ccaff030SJeremy L Thompson   case NS_ADVECTION:
1123ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1124ccaff030SJeremy L Thompson     if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxAdvection2d,
1125ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1126ccaff030SJeremy L Thompson     if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxAdvection2d,
1127ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1128ccaff030SJeremy L Thompson   }
1129ccaff030SJeremy L Thompson 
1130ccaff030SJeremy L Thompson   // Set up PETSc context
1131ccaff030SJeremy L Thompson   // Set up units structure
1132ccaff030SJeremy L Thompson   units->meter = meter;
1133ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1134ccaff030SJeremy L Thompson   units->second = second;
1135ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1136ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1137ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1138ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1139ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1140ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1141ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1142ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1143ccaff030SJeremy L Thompson 
1144ccaff030SJeremy L Thompson   // Set up user structure
1145ccaff030SJeremy L Thompson   user->comm = comm;
1146ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1147ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1148ccaff030SJeremy L Thompson   user->units = units;
1149ccaff030SJeremy L Thompson   user->dm = dm;
1150ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1151ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1152ccaff030SJeremy L Thompson   user->ceed = ceed;
1153ccaff030SJeremy L Thompson 
1154ccaff030SJeremy L Thompson   // Calculate qdata and ICs
1155ccaff030SJeremy L Thompson   // Set up state global and local vectors
1156ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1157ccaff030SJeremy L Thompson 
1158ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1159ccaff030SJeremy L Thompson 
1160ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1161ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
1162ccaff030SJeremy L Thompson   CeedOperatorApply(op_setup, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
1163ccaff030SJeremy L Thompson   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1164ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1165ccaff030SJeremy L Thompson 
116627dd567eSvaleriabarra   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
1167ccaff030SJeremy L Thompson                              &ctxSetup, 0.0);
1168ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1169ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1170ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1171ccaff030SJeremy L Thompson     // slower execution.
1172ccaff030SJeremy L Thompson     Vec Qbc;
1173ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1174ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1175ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1176ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1177ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1178ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1179ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
1180380e0a58Svaleriabarra                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
1181380e0a58Svaleriabarra     CHKERRQ(ierr);
1182ccaff030SJeremy L Thompson   }
1183ccaff030SJeremy L Thompson 
1184ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1185ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1186ccaff030SJeremy L Thompson   // Gather initial Q values
1187ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1188ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1189ccaff030SJeremy L Thompson     PetscViewer viewer;
1190ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1191ccaff030SJeremy L Thompson     // Read input
1192ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1193ccaff030SJeremy L Thompson                          user->outputfolder);
1194ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1195ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1196ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1197ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1198ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1199ccaff030SJeremy L Thompson   }
1200ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1201ccaff030SJeremy L Thompson 
1202ccaff030SJeremy L Thompson   // Create and setup TS
1203ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1204ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1205ccaff030SJeremy L Thompson   if (implicit) {
1206ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1207ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1208ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1209ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1210ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1211ccaff030SJeremy L Thompson     }
1212ccaff030SJeremy L Thompson   } else {
1213ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1214ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1215ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1216ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1217ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1218ccaff030SJeremy L Thompson   }
1219ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1220ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1221ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
122284f4e65dSvaleriabarra   if (test) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1223ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1224ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1225ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1226ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1227ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1228ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
1229ccaff030SJeremy L Thompson     if (!test) {
1230ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1231ccaff030SJeremy L Thompson     }
1232ccaff030SJeremy L Thompson   } else { // continue from time of last output
1233ccaff030SJeremy L Thompson     PetscReal time;
1234ccaff030SJeremy L Thompson     PetscInt count;
1235ccaff030SJeremy L Thompson     PetscViewer viewer;
1236ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1237ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1238ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1239ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1240ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1241ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1242ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1243ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1244ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1245ccaff030SJeremy L Thompson   }
1246ccaff030SJeremy L Thompson   if (!test) {
1247ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1248ccaff030SJeremy L Thompson   }
1249ccaff030SJeremy L Thompson 
1250ccaff030SJeremy L Thompson   // Solve
1251ccaff030SJeremy L Thompson   start = MPI_Wtime();
1252ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1253ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1254ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1255ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1256ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1257ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
1258ccaff030SJeremy L Thompson   if (!test) {
1259ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1260380e0a58Svaleriabarra                        "Time taken for solution (sec): %g\n",
1261ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1262ccaff030SJeremy L Thompson   }
1263ccaff030SJeremy L Thompson 
1264ccaff030SJeremy L Thompson   // Get error
1265ccaff030SJeremy L Thompson   if (problem->non_zero_time && !test) {
1266ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1267ccaff030SJeremy L Thompson     PetscReal norm;
1268ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1269ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1270ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1271ccaff030SJeremy L Thompson 
127227dd567eSvaleriabarra     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
1273ccaff030SJeremy L Thompson                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1274ccaff030SJeremy L Thompson 
1275ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1276ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1277ccaff030SJeremy L Thompson     CeedVectorDestroy(&q0ceed);
1278ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1279ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1280ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
1281357222b3Svaleriabarra     // Clean up vectors
1282357222b3Svaleriabarra     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1283357222b3Svaleriabarra     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1284ccaff030SJeremy L Thompson   }
1285ccaff030SJeremy L Thompson 
1286ccaff030SJeremy L Thompson   // Output Statistics
1287ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
1288ccaff030SJeremy L Thompson   if (!test) {
1289ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1290ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1291ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1292ccaff030SJeremy L Thompson   }
1293d5424b73Svaleriabarra   // Output numerical values from command line
1294d5424b73Svaleriabarra   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
1295ccaff030SJeremy L Thompson 
1296*9cf88b28Svaleriabarra   if (test) { // compare reference solution values with current run
1297*9cf88b28Svaleriabarra     PetscViewer viewer;
1298*9cf88b28Svaleriabarra     char filepath[PETSC_MAX_PATH_LEN];
1299*9cf88b28Svaleriabarra     PetscReal testtol;
1300*9cf88b28Svaleriabarra     // Create reference solution file name and assign relative error tolerance
1301*9cf88b28Svaleriabarra     if (!implicit) { // Explicit test run
1302*9cf88b28Svaleriabarra       testtol = 1E-5;
1303*9cf88b28Svaleriabarra       ierr = PetscStrncpy(filepath,
1304*9cf88b28Svaleriabarra                           "examples/fluids/tests-output/fluids-navierstokes-explicit.bin",
1305*9cf88b28Svaleriabarra                           sizeof filepath);
1306*9cf88b28Svaleriabarra       CHKERRQ(ierr);
1307*9cf88b28Svaleriabarra     }
1308*9cf88b28Svaleriabarra     else { // Implicit test runs
1309*9cf88b28Svaleriabarra       testtol = 5E-4;
1310*9cf88b28Svaleriabarra       if (stab != STAB_SUPG) {
1311*9cf88b28Svaleriabarra         ierr = PetscStrncpy(filepath,
1312*9cf88b28Svaleriabarra                             "examples/fluids/tests-output/fluids-navierstokes-implicit-nostab.bin",
1313*9cf88b28Svaleriabarra                             sizeof filepath);
1314*9cf88b28Svaleriabarra         CHKERRQ(ierr);
1315*9cf88b28Svaleriabarra       }
1316*9cf88b28Svaleriabarra       else {
1317*9cf88b28Svaleriabarra         ierr = PetscStrncpy(filepath,
1318*9cf88b28Svaleriabarra                             "examples/fluids/tests-output/fluids-navierstokes-implicit-supgstab.bin",
1319*9cf88b28Svaleriabarra                             sizeof filepath);
1320*9cf88b28Svaleriabarra         CHKERRQ(ierr);
1321*9cf88b28Svaleriabarra       }
1322*9cf88b28Svaleriabarra     }
1323*9cf88b28Svaleriabarra     // Read reference file
1324*9cf88b28Svaleriabarra     Vec Qref;
1325*9cf88b28Svaleriabarra     PetscReal error, Qrefnorm;
1326*9cf88b28Svaleriabarra     ierr = VecDuplicate(Q, &Qref);
1327*9cf88b28Svaleriabarra     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1328*9cf88b28Svaleriabarra     CHKERRQ(ierr);
1329*9cf88b28Svaleriabarra     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
1330*9cf88b28Svaleriabarra     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1331*9cf88b28Svaleriabarra 
1332*9cf88b28Svaleriabarra     // Compute error with respect to reference solution
1333*9cf88b28Svaleriabarra     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
1334*9cf88b28Svaleriabarra     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
1335*9cf88b28Svaleriabarra     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
1336*9cf88b28Svaleriabarra     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
1337*9cf88b28Svaleriabarra     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
1338*9cf88b28Svaleriabarra     // Check error
1339*9cf88b28Svaleriabarra     if (error > testtol) {
1340*9cf88b28Svaleriabarra       ierr = PetscPrintf(PETSC_COMM_WORLD,
1341*9cf88b28Svaleriabarra                          "Test failed with error norm %g\n",
1342*9cf88b28Svaleriabarra                          (double)error); CHKERRQ(ierr);
1343*9cf88b28Svaleriabarra     }
1344*9cf88b28Svaleriabarra   }
1345*9cf88b28Svaleriabarra 
1346ccaff030SJeremy L Thompson   // Clean up libCEED
1347ccaff030SJeremy L Thompson   CeedVectorDestroy(&qdata);
1348ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1349ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1350ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1351ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
1352ccaff030SJeremy L Thompson   CeedBasisDestroy(&basisq);
1353ccaff030SJeremy L Thompson   CeedBasisDestroy(&basisx);
1354ccaff030SJeremy L Thompson   CeedBasisDestroy(&basisxc);
1355ccaff030SJeremy L Thompson   CeedElemRestrictionDestroy(&restrictq);
1356ccaff030SJeremy L Thompson   CeedElemRestrictionDestroy(&restrictx);
1357ccaff030SJeremy L Thompson   CeedElemRestrictionDestroy(&restrictqdi);
1358ccaff030SJeremy L Thompson   CeedElemRestrictionDestroy(&restrictxcoord);
1359ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_setup);
1360ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1361ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_rhs);
1362ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ifunction);
1363ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_setup);
1364ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
1365ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1366ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1367ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
1368ccaff030SJeremy L Thompson 
1369ccaff030SJeremy L Thompson   // Clean up PETSc
1370ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1371ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1372ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1373ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1374ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1375ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1376ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1377ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1378ccaff030SJeremy L Thompson   return PetscFinalize();
1379ccaff030SJeremy L Thompson }
1380ccaff030SJeremy L Thompson 
1381