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