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