xref: /petsc/src/ts/tutorials/ex14.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static const char help[] = "Toy hydrostatic ice flow with multigrid in 3D.\n\
2c4762a1bSJed Brown \n\
3c4762a1bSJed Brown Solves the hydrostatic (aka Blatter/Pattyn/First Order) equations for ice sheet flow\n\
4c4762a1bSJed Brown using multigrid.  The ice uses a power-law rheology with \"Glen\" exponent 3 (corresponds\n\
5c4762a1bSJed Brown to p=4/3 in a p-Laplacian).  The focus is on ISMIP-HOM experiments which assume periodic\n\
6c4762a1bSJed Brown boundary conditions in the x- and y-directions.\n\
7c4762a1bSJed Brown \n\
8c4762a1bSJed Brown Equations are rescaled so that the domain size and solution are O(1), details of this scaling\n\
9c4762a1bSJed Brown can be controlled by the options -units_meter, -units_second, and -units_kilogram.\n\
10c4762a1bSJed Brown \n\
11c4762a1bSJed Brown A VTK StructuredGrid output file can be written using the option -o filename.vts\n\
12c4762a1bSJed Brown \n\n";
13c4762a1bSJed Brown 
14c4762a1bSJed Brown /*
15c4762a1bSJed Brown The equations for horizontal velocity (u,v) are
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   - [eta (4 u_x + 2 v_y)]_x - [eta (u_y + v_x)]_y - [eta u_z]_z + rho g s_x = 0
18c4762a1bSJed Brown   - [eta (4 v_y + 2 u_x)]_y - [eta (u_y + v_x)]_x - [eta v_z]_z + rho g s_y = 0
19c4762a1bSJed Brown 
20c4762a1bSJed Brown where
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   eta = B/2 (epsilon + gamma)^((p-2)/2)
23c4762a1bSJed Brown 
24c4762a1bSJed Brown is the nonlinear effective viscosity with regularization epsilon and hardness parameter B,
25c4762a1bSJed Brown written in terms of the second invariant
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   gamma = u_x^2 + v_y^2 + u_x v_y + (1/4) (u_y + v_x)^2 + (1/4) u_z^2 + (1/4) v_z^2
28c4762a1bSJed Brown 
29c4762a1bSJed Brown The surface boundary conditions are the natural conditions.  The basal boundary conditions
30c4762a1bSJed Brown are either no-slip, or Navier (linear) slip with spatially variant friction coefficient beta^2.
31c4762a1bSJed Brown 
32c4762a1bSJed Brown In the code, the equations for (u,v) are multiplied through by 1/(rho g) so that residuals are O(1).
33c4762a1bSJed Brown 
34c4762a1bSJed Brown The discretization is Q1 finite elements, managed by a DMDA.  The grid is never distorted in the
35c4762a1bSJed Brown map (x,y) plane, but the bed and surface may be bumpy.  This is handled as usual in FEM, through
36c4762a1bSJed Brown the Jacobian of the coordinate transformation from a reference element to the physical element.
37c4762a1bSJed Brown 
38c4762a1bSJed Brown Since ice-flow is tightly coupled in the z-direction (within columns), the DMDA is managed
39c4762a1bSJed Brown specially so that columns are never distributed, and are always contiguous in memory.
40c4762a1bSJed Brown This amounts to reversing the meaning of X,Y,Z compared to the DMDA's internal interpretation,
41c4762a1bSJed Brown and then indexing as vec[i][j][k].  The exotic coarse spaces require 2D DMDAs which are made to
42c4762a1bSJed Brown use compatible domain decomposition relative to the 3D DMDAs.
43c4762a1bSJed Brown 
44c4762a1bSJed Brown */
45c4762a1bSJed Brown 
46c4762a1bSJed Brown #include <petscts.h>
47c4762a1bSJed Brown #include <petscdm.h>
48c4762a1bSJed Brown #include <petscdmda.h>
49c4762a1bSJed Brown #include <petscdmcomposite.h>
50c4762a1bSJed Brown #include <ctype.h> /* toupper() */
51c4762a1bSJed Brown #include <petsc/private/petscimpl.h>
52c4762a1bSJed Brown 
53c4762a1bSJed Brown #if defined __SSE2__
54c4762a1bSJed Brown #include <emmintrin.h>
55c4762a1bSJed Brown #endif
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /* The SSE2 kernels are only for PetscScalar=double on architectures that support it */
589371c9d4SSatish Balay #define USE_SSE2_KERNELS (!defined NO_SSE2 && !defined PETSC_USE_COMPLEX && !defined PETSC_USE_REAL_SINGLE && defined __SSE2__)
59c4762a1bSJed Brown 
60c4762a1bSJed Brown #if !defined __STDC_VERSION__ || __STDC_VERSION__ < 199901L
61c4762a1bSJed Brown #if defined  __cplusplus /* C++ restrict is nonstandard and compilers have inconsistent rules about where it can be used */
62c4762a1bSJed Brown #define restrict
63c4762a1bSJed Brown #else
64c4762a1bSJed Brown #define restrict PETSC_RESTRICT
65c4762a1bSJed Brown #endif
66c4762a1bSJed Brown #endif
67c4762a1bSJed Brown 
68c4762a1bSJed Brown static PetscClassId THI_CLASSID;
69c4762a1bSJed Brown 
709371c9d4SSatish Balay typedef enum {
719371c9d4SSatish Balay   QUAD_GAUSS,
729371c9d4SSatish Balay   QUAD_LOBATTO
739371c9d4SSatish Balay } QuadratureType;
74c4762a1bSJed Brown static const char     *QuadratureTypes[] = {"gauss", "lobatto", "QuadratureType", "QUAD_", 0};
75c4762a1bSJed Brown static const PetscReal HexQWeights[8]    = {1, 1, 1, 1, 1, 1, 1, 1};
76c4762a1bSJed Brown static const PetscReal HexQNodes[]       = {-0.57735026918962573, 0.57735026918962573};
77c4762a1bSJed Brown #define G 0.57735026918962573
78c4762a1bSJed Brown #define H (0.5 * (1. + G))
79c4762a1bSJed Brown #define L (0.5 * (1. - G))
80c4762a1bSJed Brown #define M (-0.5)
81c4762a1bSJed Brown #define P (0.5)
82c4762a1bSJed Brown /* Special quadrature: Lobatto in horizontal, Gauss in vertical */
839371c9d4SSatish Balay static const PetscReal HexQInterp_Lobatto[8][8] = {
849371c9d4SSatish Balay   {H, 0, 0, 0, L, 0, 0, 0},
85c4762a1bSJed Brown   {0, H, 0, 0, 0, L, 0, 0},
86c4762a1bSJed Brown   {0, 0, H, 0, 0, 0, L, 0},
87c4762a1bSJed Brown   {0, 0, 0, H, 0, 0, 0, L},
88c4762a1bSJed Brown   {L, 0, 0, 0, H, 0, 0, 0},
89c4762a1bSJed Brown   {0, L, 0, 0, 0, H, 0, 0},
90c4762a1bSJed Brown   {0, 0, L, 0, 0, 0, H, 0},
919371c9d4SSatish Balay   {0, 0, 0, L, 0, 0, 0, H}
929371c9d4SSatish Balay };
93c4762a1bSJed Brown static const PetscReal HexQDeriv_Lobatto[8][8][3] = {
94c4762a1bSJed Brown   {{M * H, M *H, M}, {P * H, 0, 0},    {0, 0, 0},        {0, P *H, 0},     {M * L, M *L, P}, {P * L, 0, 0},    {0, 0, 0},        {0, P *L, 0}    },
95c4762a1bSJed Brown   {{M * H, 0, 0},    {P * H, M *H, M}, {0, P *H, 0},     {0, 0, 0},        {M * L, 0, 0},    {P * L, M *L, P}, {0, P *L, 0},     {0, 0, 0}       },
96c4762a1bSJed Brown   {{0, 0, 0},        {0, M *H, 0},     {P * H, P *H, M}, {M * H, 0, 0},    {0, 0, 0},        {0, M *L, 0},     {P * L, P *L, P}, {M * L, 0, 0}   },
97c4762a1bSJed Brown   {{0, M *H, 0},     {0, 0, 0},        {P * H, 0, 0},    {M * H, P *H, M}, {0, M *L, 0},     {0, 0, 0},        {P * L, 0, 0},    {M * L, P *L, P}},
98c4762a1bSJed Brown   {{M * L, M *L, M}, {P * L, 0, 0},    {0, 0, 0},        {0, P *L, 0},     {M * H, M *H, P}, {P * H, 0, 0},    {0, 0, 0},        {0, P *H, 0}    },
99c4762a1bSJed Brown   {{M * L, 0, 0},    {P * L, M *L, M}, {0, P *L, 0},     {0, 0, 0},        {M * H, 0, 0},    {P * H, M *H, P}, {0, P *H, 0},     {0, 0, 0}       },
100c4762a1bSJed Brown   {{0, 0, 0},        {0, M *L, 0},     {P * L, P *L, M}, {M * L, 0, 0},    {0, 0, 0},        {0, M *H, 0},     {P * H, P *H, P}, {M * H, 0, 0}   },
1019371c9d4SSatish Balay   {{0, M *L, 0},     {0, 0, 0},        {P * L, 0, 0},    {M * L, P *L, M}, {0, M *H, 0},     {0, 0, 0},        {P * H, 0, 0},    {M * H, P *H, P}}
1029371c9d4SSatish Balay };
103c4762a1bSJed Brown /* Stanndard Gauss */
1049371c9d4SSatish Balay static const PetscReal HexQInterp_Gauss[8][8] = {
1059371c9d4SSatish Balay   {H * H * H, L *H *H, L *L *H, H *L *H, H *H *L, L *H *L, L *L *L, H *L *L},
106c4762a1bSJed Brown   {L * H * H, H *H *H, H *L *H, L *L *H, L *H *L, H *H *L, H *L *L, L *L *L},
107c4762a1bSJed Brown   {L * L * H, H *L *H, H *H *H, L *H *H, L *L *L, H *L *L, H *H *L, L *H *L},
108c4762a1bSJed Brown   {H * L * H, L *L *H, L *H *H, H *H *H, H *L *L, L *L *L, L *H *L, H *H *L},
109c4762a1bSJed Brown   {H * H * L, L *H *L, L *L *L, H *L *L, H *H *H, L *H *H, L *L *H, H *L *H},
110c4762a1bSJed Brown   {L * H * L, H *H *L, H *L *L, L *L *L, L *H *H, H *H *H, H *L *H, L *L *H},
111c4762a1bSJed Brown   {L * L * L, H *L *L, H *H *L, L *H *L, L *L *H, H *L *H, H *H *H, L *H *H},
1129371c9d4SSatish Balay   {H * L * L, L *L *L, L *H *L, H *H *L, H *L *H, L *L *H, L *H *H, H *H *H}
1139371c9d4SSatish Balay };
114c4762a1bSJed Brown static const PetscReal HexQDeriv_Gauss[8][8][3] = {
115c4762a1bSJed Brown   {{M * H * H, H *M *H, H *H *M}, {P * H * H, L *M *H, L *H *M}, {P * L * H, L *P *H, L *L *M}, {M * L * H, H *P *H, H *L *M}, {M * H * L, H *M *L, H *H *P}, {P * H * L, L *M *L, L *H *P}, {P * L * L, L *P *L, L *L *P}, {M * L * L, H *P *L, H *L *P}},
116c4762a1bSJed Brown   {{M * H * H, L *M *H, L *H *M}, {P * H * H, H *M *H, H *H *M}, {P * L * H, H *P *H, H *L *M}, {M * L * H, L *P *H, L *L *M}, {M * H * L, L *M *L, L *H *P}, {P * H * L, H *M *L, H *H *P}, {P * L * L, H *P *L, H *L *P}, {M * L * L, L *P *L, L *L *P}},
117c4762a1bSJed Brown   {{M * L * H, L *M *H, L *L *M}, {P * L * H, H *M *H, H *L *M}, {P * H * H, H *P *H, H *H *M}, {M * H * H, L *P *H, L *H *M}, {M * L * L, L *M *L, L *L *P}, {P * L * L, H *M *L, H *L *P}, {P * H * L, H *P *L, H *H *P}, {M * H * L, L *P *L, L *H *P}},
118c4762a1bSJed Brown   {{M * L * H, H *M *H, H *L *M}, {P * L * H, L *M *H, L *L *M}, {P * H * H, L *P *H, L *H *M}, {M * H * H, H *P *H, H *H *M}, {M * L * L, H *M *L, H *L *P}, {P * L * L, L *M *L, L *L *P}, {P * H * L, L *P *L, L *H *P}, {M * H * L, H *P *L, H *H *P}},
119c4762a1bSJed Brown   {{M * H * L, H *M *L, H *H *M}, {P * H * L, L *M *L, L *H *M}, {P * L * L, L *P *L, L *L *M}, {M * L * L, H *P *L, H *L *M}, {M * H * H, H *M *H, H *H *P}, {P * H * H, L *M *H, L *H *P}, {P * L * H, L *P *H, L *L *P}, {M * L * H, H *P *H, H *L *P}},
120c4762a1bSJed Brown   {{M * H * L, L *M *L, L *H *M}, {P * H * L, H *M *L, H *H *M}, {P * L * L, H *P *L, H *L *M}, {M * L * L, L *P *L, L *L *M}, {M * H * H, L *M *H, L *H *P}, {P * H * H, H *M *H, H *H *P}, {P * L * H, H *P *H, H *L *P}, {M * L * H, L *P *H, L *L *P}},
121c4762a1bSJed Brown   {{M * L * L, L *M *L, L *L *M}, {P * L * L, H *M *L, H *L *M}, {P * H * L, H *P *L, H *H *M}, {M * H * L, L *P *L, L *H *M}, {M * L * H, L *M *H, L *L *P}, {P * L * H, H *M *H, H *L *P}, {P * H * H, H *P *H, H *H *P}, {M * H * H, L *P *H, L *H *P}},
1229371c9d4SSatish Balay   {{M * L * L, H *M *L, H *L *M}, {P * L * L, L *M *L, L *L *M}, {P * H * L, L *P *L, L *H *M}, {M * H * L, H *P *L, H *H *M}, {M * L * H, H *M *H, H *L *P}, {P * L * H, L *M *H, L *L *P}, {P * H * H, L *P *H, L *H *P}, {M * H * H, H *P *H, H *H *P}}
1239371c9d4SSatish Balay };
124c4762a1bSJed Brown static const PetscReal (*HexQInterp)[8], (*HexQDeriv)[8][3];
125c4762a1bSJed Brown /* Standard 2x2 Gauss quadrature for the bottom layer. */
1269371c9d4SSatish Balay static const PetscReal QuadQInterp[4][4] = {
1279371c9d4SSatish Balay   {H * H, L *H, L *L, H *L},
128c4762a1bSJed Brown   {L * H, H *H, H *L, L *L},
129c4762a1bSJed Brown   {L * L, H *L, H *H, L *H},
1309371c9d4SSatish Balay   {H * L, L *L, L *H, H *H}
1319371c9d4SSatish Balay };
132c4762a1bSJed Brown static const PetscReal QuadQDeriv[4][4][2] = {
133c4762a1bSJed Brown   {{M * H, M *H}, {P * H, M *L}, {P * L, P *L}, {M * L, P *H}},
134c4762a1bSJed Brown   {{M * H, M *L}, {P * H, M *H}, {P * L, P *H}, {M * L, P *L}},
135c4762a1bSJed Brown   {{M * L, M *L}, {P * L, M *H}, {P * H, P *H}, {M * H, P *L}},
1369371c9d4SSatish Balay   {{M * L, M *H}, {P * L, M *L}, {P * H, P *L}, {M * H, P *H}}
1379371c9d4SSatish Balay };
138c4762a1bSJed Brown #undef G
139c4762a1bSJed Brown #undef H
140c4762a1bSJed Brown #undef L
141c4762a1bSJed Brown #undef M
142c4762a1bSJed Brown #undef P
143c4762a1bSJed Brown 
1449371c9d4SSatish Balay #define HexExtract(x, i, j, k, n) \
1459371c9d4SSatish Balay   do { \
146c4762a1bSJed Brown     (n)[0] = (x)[i][j][k]; \
147c4762a1bSJed Brown     (n)[1] = (x)[i + 1][j][k]; \
148c4762a1bSJed Brown     (n)[2] = (x)[i + 1][j + 1][k]; \
149c4762a1bSJed Brown     (n)[3] = (x)[i][j + 1][k]; \
150c4762a1bSJed Brown     (n)[4] = (x)[i][j][k + 1]; \
151c4762a1bSJed Brown     (n)[5] = (x)[i + 1][j][k + 1]; \
152c4762a1bSJed Brown     (n)[6] = (x)[i + 1][j + 1][k + 1]; \
153c4762a1bSJed Brown     (n)[7] = (x)[i][j + 1][k + 1]; \
154c4762a1bSJed Brown   } while (0)
155c4762a1bSJed Brown 
1569371c9d4SSatish Balay #define HexExtractRef(x, i, j, k, n) \
1579371c9d4SSatish Balay   do { \
158c4762a1bSJed Brown     (n)[0] = &(x)[i][j][k]; \
159c4762a1bSJed Brown     (n)[1] = &(x)[i + 1][j][k]; \
160c4762a1bSJed Brown     (n)[2] = &(x)[i + 1][j + 1][k]; \
161c4762a1bSJed Brown     (n)[3] = &(x)[i][j + 1][k]; \
162c4762a1bSJed Brown     (n)[4] = &(x)[i][j][k + 1]; \
163c4762a1bSJed Brown     (n)[5] = &(x)[i + 1][j][k + 1]; \
164c4762a1bSJed Brown     (n)[6] = &(x)[i + 1][j + 1][k + 1]; \
165c4762a1bSJed Brown     (n)[7] = &(x)[i][j + 1][k + 1]; \
166c4762a1bSJed Brown   } while (0)
167c4762a1bSJed Brown 
1689371c9d4SSatish Balay #define QuadExtract(x, i, j, n) \
1699371c9d4SSatish Balay   do { \
170c4762a1bSJed Brown     (n)[0] = (x)[i][j]; \
171c4762a1bSJed Brown     (n)[1] = (x)[i + 1][j]; \
172c4762a1bSJed Brown     (n)[2] = (x)[i + 1][j + 1]; \
173c4762a1bSJed Brown     (n)[3] = (x)[i][j + 1]; \
174c4762a1bSJed Brown   } while (0)
175c4762a1bSJed Brown 
1769371c9d4SSatish Balay static PetscScalar Sqr(PetscScalar a) {
1779371c9d4SSatish Balay   return a * a;
1789371c9d4SSatish Balay }
179c4762a1bSJed Brown 
1809371c9d4SSatish Balay static void HexGrad(const PetscReal dphi[][3], const PetscReal zn[], PetscReal dz[]) {
181c4762a1bSJed Brown   PetscInt i;
182c4762a1bSJed Brown   dz[0] = dz[1] = dz[2] = 0;
183c4762a1bSJed Brown   for (i = 0; i < 8; i++) {
184c4762a1bSJed Brown     dz[0] += dphi[i][0] * zn[i];
185c4762a1bSJed Brown     dz[1] += dphi[i][1] * zn[i];
186c4762a1bSJed Brown     dz[2] += dphi[i][2] * zn[i];
187c4762a1bSJed Brown   }
188c4762a1bSJed Brown }
189c4762a1bSJed Brown 
1909371c9d4SSatish Balay static void HexComputeGeometry(PetscInt q, PetscReal hx, PetscReal hy, const PetscReal dz[restrict], PetscReal phi[restrict], PetscReal dphi[restrict][3], PetscReal *restrict jw) {
1919371c9d4SSatish Balay   const PetscReal jac[3][3] =
192c4762a1bSJed Brown     {
1939371c9d4SSatish Balay       {hx / 2, 0,      0    },
1949371c9d4SSatish Balay       {0,      hy / 2, 0    },
1959371c9d4SSatish Balay       {dz[0],  dz[1],  dz[2]}
1969371c9d4SSatish Balay   },
1979371c9d4SSatish Balay                   ijac[3][3] = {{1 / jac[0][0], 0, 0}, {0, 1 / jac[1][1], 0}, {-jac[2][0] / (jac[0][0] * jac[2][2]), -jac[2][1] / (jac[1][1] * jac[2][2]), 1 / jac[2][2]}}, jdet = jac[0][0] * jac[1][1] * jac[2][2];
198c4762a1bSJed Brown   PetscInt i;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   for (i = 0; i < 8; i++) {
201c4762a1bSJed Brown     const PetscReal *dphir = HexQDeriv[q][i];
202c4762a1bSJed Brown     phi[i]                 = HexQInterp[q][i];
203c4762a1bSJed Brown     dphi[i][0]             = dphir[0] * ijac[0][0] + dphir[1] * ijac[1][0] + dphir[2] * ijac[2][0];
204c4762a1bSJed Brown     dphi[i][1]             = dphir[0] * ijac[0][1] + dphir[1] * ijac[1][1] + dphir[2] * ijac[2][1];
205c4762a1bSJed Brown     dphi[i][2]             = dphir[0] * ijac[0][2] + dphir[1] * ijac[1][2] + dphir[2] * ijac[2][2];
206c4762a1bSJed Brown   }
207c4762a1bSJed Brown   *jw = 1.0 * jdet;
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown typedef struct _p_THI   *THI;
211c4762a1bSJed Brown typedef struct _n_Units *Units;
212c4762a1bSJed Brown 
213c4762a1bSJed Brown typedef struct {
214c4762a1bSJed Brown   PetscScalar u, v;
215c4762a1bSJed Brown } Node;
216c4762a1bSJed Brown 
217c4762a1bSJed Brown typedef struct {
218c4762a1bSJed Brown   PetscScalar b;     /* bed */
219c4762a1bSJed Brown   PetscScalar h;     /* thickness */
220c4762a1bSJed Brown   PetscScalar beta2; /* friction */
221c4762a1bSJed Brown } PrmNode;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown #define FieldSize(ntype)             ((PetscInt)(sizeof(ntype) / sizeof(PetscScalar)))
224c4762a1bSJed Brown #define FieldOffset(ntype, member)   ((PetscInt)(offsetof(ntype, member) / sizeof(PetscScalar)))
225c4762a1bSJed Brown #define FieldIndex(ntype, i, member) ((PetscInt)((i)*FieldSize(ntype) + FieldOffset(ntype, member)))
226c4762a1bSJed Brown #define NODE_SIZE                    FieldSize(Node)
227c4762a1bSJed Brown #define PRMNODE_SIZE                 FieldSize(PrmNode)
228c4762a1bSJed Brown 
229c4762a1bSJed Brown typedef struct {
230c4762a1bSJed Brown   PetscReal min, max, cmin, cmax;
231c4762a1bSJed Brown } PRange;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown struct _p_THI {
234c4762a1bSJed Brown   PETSCHEADER(int);
235c4762a1bSJed Brown   void (*initialize)(THI, PetscReal x, PetscReal y, PrmNode *p);
236c4762a1bSJed Brown   PetscInt  nlevels;
237c4762a1bSJed Brown   PetscInt  zlevels;
238c4762a1bSJed Brown   PetscReal Lx, Ly, Lz; /* Model domain */
239c4762a1bSJed Brown   PetscReal alpha;      /* Bed angle */
240c4762a1bSJed Brown   Units     units;
241c4762a1bSJed Brown   PetscReal dirichlet_scale;
242c4762a1bSJed Brown   PetscReal ssa_friction_scale;
243c4762a1bSJed Brown   PetscReal inertia;
244c4762a1bSJed Brown   PRange    eta;
245c4762a1bSJed Brown   PRange    beta2;
246c4762a1bSJed Brown   struct {
247c4762a1bSJed Brown     PetscReal Bd2, eps, exponent, glen_n;
248c4762a1bSJed Brown   } viscosity;
249c4762a1bSJed Brown   struct {
250c4762a1bSJed Brown     PetscReal irefgam, eps2, exponent;
251c4762a1bSJed Brown   } friction;
252c4762a1bSJed Brown   struct {
253c4762a1bSJed Brown     PetscReal rate, exponent, refvel;
254c4762a1bSJed Brown   } erosion;
255c4762a1bSJed Brown   PetscReal rhog;
256c4762a1bSJed Brown   PetscBool no_slip;
257c4762a1bSJed Brown   PetscBool verbose;
258c4762a1bSJed Brown   char     *mattype;
259c4762a1bSJed Brown   char     *monitor_basename;
260c4762a1bSJed Brown   PetscInt  monitor_interval;
261c4762a1bSJed Brown };
262c4762a1bSJed Brown 
263c4762a1bSJed Brown struct _n_Units {
264c4762a1bSJed Brown   /* fundamental */
265c4762a1bSJed Brown   PetscReal meter;
266c4762a1bSJed Brown   PetscReal kilogram;
267c4762a1bSJed Brown   PetscReal second;
268c4762a1bSJed Brown   /* derived */
269c4762a1bSJed Brown   PetscReal Pascal;
270c4762a1bSJed Brown   PetscReal year;
271c4762a1bSJed Brown };
272c4762a1bSJed Brown 
2739371c9d4SSatish Balay static void PrmHexGetZ(const PrmNode pn[], PetscInt k, PetscInt zm, PetscReal zn[]) {
2749371c9d4SSatish Balay   const PetscScalar zm1 = zm - 1, znl[8] = {pn[0].b + pn[0].h * (PetscScalar)k / zm1,       pn[1].b + pn[1].h * (PetscScalar)k / zm1,       pn[2].b + pn[2].h * (PetscScalar)k / zm1,       pn[3].b + pn[3].h * (PetscScalar)k / zm1,
2759371c9d4SSatish Balay                                             pn[0].b + pn[0].h * (PetscScalar)(k + 1) / zm1, pn[1].b + pn[1].h * (PetscScalar)(k + 1) / zm1, pn[2].b + pn[2].h * (PetscScalar)(k + 1) / zm1, pn[3].b + pn[3].h * (PetscScalar)(k + 1) / zm1};
276c4762a1bSJed Brown   PetscInt          i;
277c4762a1bSJed Brown   for (i = 0; i < 8; i++) zn[i] = PetscRealPart(znl[i]);
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280c4762a1bSJed Brown /* Compute a gradient of all the 2D fields at four quadrature points.  Output for [quadrature_point][direction].field_name */
2819371c9d4SSatish Balay static PetscErrorCode QuadComputeGrad4(const PetscReal dphi[][4][2], PetscReal hx, PetscReal hy, const PrmNode pn[4], PrmNode dp[4][2]) {
282c4762a1bSJed Brown   PetscInt q, i, f;
283c4762a1bSJed Brown   const PetscScalar(*restrict pg)[PRMNODE_SIZE] = (const PetscScalar(*)[PRMNODE_SIZE])pn; /* Get generic array pointers to the node */
284c4762a1bSJed Brown   PetscScalar(*restrict dpg)[2][PRMNODE_SIZE]   = (PetscScalar(*)[2][PRMNODE_SIZE])dp;
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   PetscFunctionBeginUser;
2879566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(dpg, 4));
288c4762a1bSJed Brown   for (q = 0; q < 4; q++) {
289c4762a1bSJed Brown     for (i = 0; i < 4; i++) {
290c4762a1bSJed Brown       for (f = 0; f < PRMNODE_SIZE; f++) {
291c4762a1bSJed Brown         dpg[q][0][f] += dphi[q][i][0] / hx * pg[i][f];
292c4762a1bSJed Brown         dpg[q][1][f] += dphi[q][i][1] / hy * pg[i][f];
293c4762a1bSJed Brown       }
294c4762a1bSJed Brown     }
295c4762a1bSJed Brown   }
296c4762a1bSJed Brown   PetscFunctionReturn(0);
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
2999371c9d4SSatish Balay static inline PetscReal StaggeredMidpoint2D(PetscScalar a, PetscScalar b, PetscScalar c, PetscScalar d) {
3009371c9d4SSatish Balay   return 0.5 * PetscRealPart(0.75 * a + 0.75 * b + 0.25 * c + 0.25 * d);
3019371c9d4SSatish Balay }
3029371c9d4SSatish Balay static inline PetscReal UpwindFlux1D(PetscReal u, PetscReal hL, PetscReal hR) {
3039371c9d4SSatish Balay   return (u > 0) ? hL * u : hR * u;
3049371c9d4SSatish Balay }
305c4762a1bSJed Brown 
3069371c9d4SSatish Balay #define UpwindFluxXW(x3, x2, h, i, j, k, dj) \
3079371c9d4SSatish Balay   UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].u, x3[i - 1][j][k].u, x3[i - 1][j + dj][k].u, x3[i][k + dj][k].u), PetscRealPart(0.75 * x2[i - 1][j].h + 0.25 * x2[i - 1][j + dj].h), PetscRealPart(0.75 * x2[i][j].h + 0.25 * x2[i][j + dj].h))
3089371c9d4SSatish Balay #define UpwindFluxXE(x3, x2, h, i, j, k, dj) \
3099371c9d4SSatish Balay   UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].u, x3[i + 1][j][k].u, x3[i + 1][j + dj][k].u, x3[i][k + dj][k].u), PetscRealPart(0.75 * x2[i][j].h + 0.25 * x2[i][j + dj].h), PetscRealPart(0.75 * x2[i + 1][j].h + 0.25 * x2[i + 1][j + dj].h))
3109371c9d4SSatish Balay #define UpwindFluxYS(x3, x2, h, i, j, k, di) \
3119371c9d4SSatish Balay   UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].v, x3[i][j - 1][k].v, x3[i + di][j - 1][k].v, x3[i + di][j][k].v), PetscRealPart(0.75 * x2[i][j - 1].h + 0.25 * x2[i + di][j - 1].h), PetscRealPart(0.75 * x2[i][j].h + 0.25 * x2[i + di][j].h))
3129371c9d4SSatish Balay #define UpwindFluxYN(x3, x2, h, i, j, k, di) \
3139371c9d4SSatish Balay   UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].v, x3[i][j + 1][k].v, x3[i + di][j + 1][k].v, x3[i + di][j][k].v), PetscRealPart(0.75 * x2[i][j].h + 0.25 * x2[i + di][j].h), PetscRealPart(0.75 * x2[i][j + 1].h + 0.25 * x2[i + di][j + 1].h))
314c4762a1bSJed Brown 
3159371c9d4SSatish Balay static void PrmNodeGetFaceMeasure(const PrmNode **p, PetscInt i, PetscInt j, PetscScalar h[]) {
316c4762a1bSJed Brown   /* West */
317c4762a1bSJed Brown   h[0] = StaggeredMidpoint2D(p[i][j].h, p[i - 1][j].h, p[i - 1][j - 1].h, p[i][j - 1].h);
318c4762a1bSJed Brown   h[1] = StaggeredMidpoint2D(p[i][j].h, p[i - 1][j].h, p[i - 1][j + 1].h, p[i][j + 1].h);
319c4762a1bSJed Brown   /* East */
320c4762a1bSJed Brown   h[2] = StaggeredMidpoint2D(p[i][j].h, p[i + 1][j].h, p[i + 1][j + 1].h, p[i][j + 1].h);
321c4762a1bSJed Brown   h[3] = StaggeredMidpoint2D(p[i][j].h, p[i + 1][j].h, p[i + 1][j - 1].h, p[i][j - 1].h);
322c4762a1bSJed Brown   /* South */
323c4762a1bSJed Brown   h[4] = StaggeredMidpoint2D(p[i][j].h, p[i][j - 1].h, p[i + 1][j - 1].h, p[i + 1][j].h);
324c4762a1bSJed Brown   h[5] = StaggeredMidpoint2D(p[i][j].h, p[i][j - 1].h, p[i - 1][j - 1].h, p[i - 1][j].h);
325c4762a1bSJed Brown   /* North */
326c4762a1bSJed Brown   h[6] = StaggeredMidpoint2D(p[i][j].h, p[i][j + 1].h, p[i - 1][j + 1].h, p[i - 1][j].h);
327c4762a1bSJed Brown   h[7] = StaggeredMidpoint2D(p[i][j].h, p[i][j + 1].h, p[i + 1][j + 1].h, p[i + 1][j].h);
328c4762a1bSJed Brown }
329c4762a1bSJed Brown 
330c4762a1bSJed Brown /* Tests A and C are from the ISMIP-HOM paper (Pattyn et al. 2008) */
3319371c9d4SSatish Balay static void THIInitialize_HOM_A(THI thi, PetscReal x, PetscReal y, PrmNode *p) {
332c4762a1bSJed Brown   Units     units = thi->units;
333c4762a1bSJed Brown   PetscReal s     = -x * PetscSinReal(thi->alpha);
334c4762a1bSJed Brown   p->b            = s - 1000 * units->meter + 500 * units->meter * PetscSinReal(x * 2 * PETSC_PI / thi->Lx) * PetscSinReal(y * 2 * PETSC_PI / thi->Ly);
335c4762a1bSJed Brown   p->h            = s - p->b;
336c4762a1bSJed Brown   p->beta2        = -1e-10; /* This value is not used, but it should not be huge because that would change the finite difference step size  */
337c4762a1bSJed Brown }
338c4762a1bSJed Brown 
3399371c9d4SSatish Balay static void THIInitialize_HOM_C(THI thi, PetscReal x, PetscReal y, PrmNode *p) {
340c4762a1bSJed Brown   Units     units = thi->units;
341c4762a1bSJed Brown   PetscReal s     = -x * PetscSinReal(thi->alpha);
342c4762a1bSJed Brown   p->b            = s - 1000 * units->meter;
343c4762a1bSJed Brown   p->h            = s - p->b;
344c4762a1bSJed Brown   /* tau_b = beta2 v   is a stress (Pa).
345c4762a1bSJed Brown    * This is a big number in our units (it needs to balance the driving force from the surface), so we scale it by 1/rhog, just like the residual. */
346c4762a1bSJed Brown   p->beta2        = 1000 * (1 + PetscSinReal(x * 2 * PETSC_PI / thi->Lx) * PetscSinReal(y * 2 * PETSC_PI / thi->Ly)) * units->Pascal * units->year / units->meter / thi->rhog;
347c4762a1bSJed Brown }
348c4762a1bSJed Brown 
349c4762a1bSJed Brown /* These are just toys */
350c4762a1bSJed Brown 
351c4762a1bSJed Brown /* From Fred Herman */
3529371c9d4SSatish Balay static void THIInitialize_HOM_F(THI thi, PetscReal x, PetscReal y, PrmNode *p) {
353c4762a1bSJed Brown   Units     units = thi->units;
354c4762a1bSJed Brown   PetscReal s     = -x * PetscSinReal(thi->alpha);
355c4762a1bSJed Brown   p->b            = s - 1000 * units->meter + 100 * units->meter * PetscSinReal(x * 2 * PETSC_PI / thi->Lx); /* * sin(y*2*PETSC_PI/thi->Ly); */
356c4762a1bSJed Brown   p->h            = s - p->b;
357c4762a1bSJed Brown   p->h            = (1 - (atan((x - thi->Lx / 2) / 1.) + PETSC_PI / 2.) / PETSC_PI) * 500 * units->meter + 1 * units->meter;
358c4762a1bSJed Brown   s               = PetscRealPart(p->b + p->h);
359c4762a1bSJed Brown   p->beta2        = -1e-10;
360c4762a1bSJed Brown   /*  p->beta2 = 1000 * units->Pascal * units->year / units->meter; */
361c4762a1bSJed Brown }
362c4762a1bSJed Brown 
363c4762a1bSJed Brown /* Same bed as test A, free slip everywhere except for a discontinuous jump to a circular sticky region in the middle. */
3649371c9d4SSatish Balay static void THIInitialize_HOM_X(THI thi, PetscReal xx, PetscReal yy, PrmNode *p) {
365c4762a1bSJed Brown   Units     units = thi->units;
366c4762a1bSJed Brown   PetscReal x = xx * 2 * PETSC_PI / thi->Lx - PETSC_PI, y = yy * 2 * PETSC_PI / thi->Ly - PETSC_PI; /* [-pi,pi] */
367c4762a1bSJed Brown   PetscReal r = PetscSqrtReal(x * x + y * y), s = -x * PetscSinReal(thi->alpha);
368c4762a1bSJed Brown   p->b     = s - 1000 * units->meter + 500 * units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
369c4762a1bSJed Brown   p->h     = s - p->b;
370c4762a1bSJed Brown   p->beta2 = 1000 * (r < 1 ? 2 : 0) * units->Pascal * units->year / units->meter / thi->rhog;
371c4762a1bSJed Brown }
372c4762a1bSJed Brown 
373c4762a1bSJed Brown /* Like Z, but with 200 meter cliffs */
3749371c9d4SSatish Balay static void THIInitialize_HOM_Y(THI thi, PetscReal xx, PetscReal yy, PrmNode *p) {
375c4762a1bSJed Brown   Units     units = thi->units;
376c4762a1bSJed Brown   PetscReal x = xx * 2 * PETSC_PI / thi->Lx - PETSC_PI, y = yy * 2 * PETSC_PI / thi->Ly - PETSC_PI; /* [-pi,pi] */
377c4762a1bSJed Brown   PetscReal r = PetscSqrtReal(x * x + y * y), s = -x * PetscSinReal(thi->alpha);
378c4762a1bSJed Brown   p->b = s - 1000 * units->meter + 500 * units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
379c4762a1bSJed Brown   if (PetscRealPart(p->b) > -700 * units->meter) p->b += 200 * units->meter;
380c4762a1bSJed Brown   p->h     = s - p->b;
381c4762a1bSJed Brown   p->beta2 = 1000 * (1. + PetscSinReal(PetscSqrtReal(16 * r)) / PetscSqrtReal(1e-2 + 16 * r) * PetscCosReal(x * 3 / 2) * PetscCosReal(y * 3 / 2)) * units->Pascal * units->year / units->meter / thi->rhog;
382c4762a1bSJed Brown }
383c4762a1bSJed Brown 
384c4762a1bSJed Brown /* Same bed as A, smoothly varying slipperiness, similar to MATLAB's "sombrero" (uncorrelated with bathymetry) */
3859371c9d4SSatish Balay static void THIInitialize_HOM_Z(THI thi, PetscReal xx, PetscReal yy, PrmNode *p) {
386c4762a1bSJed Brown   Units     units = thi->units;
387c4762a1bSJed Brown   PetscReal x = xx * 2 * PETSC_PI / thi->Lx - PETSC_PI, y = yy * 2 * PETSC_PI / thi->Ly - PETSC_PI; /* [-pi,pi] */
388c4762a1bSJed Brown   PetscReal r = PetscSqrtReal(x * x + y * y), s = -x * PetscSinReal(thi->alpha);
389c4762a1bSJed Brown   p->b     = s - 1000 * units->meter + 500 * units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
390c4762a1bSJed Brown   p->h     = s - p->b;
391c4762a1bSJed Brown   p->beta2 = 1000 * (1. + PetscSinReal(PetscSqrtReal(16 * r)) / PetscSqrtReal(1e-2 + 16 * r) * PetscCosReal(x * 3 / 2) * PetscCosReal(y * 3 / 2)) * units->Pascal * units->year / units->meter / thi->rhog;
392c4762a1bSJed Brown }
393c4762a1bSJed Brown 
3949371c9d4SSatish Balay static void THIFriction(THI thi, PetscReal rbeta2, PetscReal gam, PetscReal *beta2, PetscReal *dbeta2) {
395c4762a1bSJed Brown   if (thi->friction.irefgam == 0) {
396c4762a1bSJed Brown     Units units           = thi->units;
397c4762a1bSJed Brown     thi->friction.irefgam = 1. / (0.5 * PetscSqr(100 * units->meter / units->year));
398c4762a1bSJed Brown     thi->friction.eps2    = 0.5 * PetscSqr(1.e-4 / thi->friction.irefgam);
399c4762a1bSJed Brown   }
400c4762a1bSJed Brown   if (thi->friction.exponent == 0) {
401c4762a1bSJed Brown     *beta2  = rbeta2;
402c4762a1bSJed Brown     *dbeta2 = 0;
403c4762a1bSJed Brown   } else {
404c4762a1bSJed Brown     *beta2  = rbeta2 * PetscPowReal(thi->friction.eps2 + gam * thi->friction.irefgam, thi->friction.exponent);
405c4762a1bSJed Brown     *dbeta2 = thi->friction.exponent * *beta2 / (thi->friction.eps2 + gam * thi->friction.irefgam) * thi->friction.irefgam;
406c4762a1bSJed Brown   }
407c4762a1bSJed Brown }
408c4762a1bSJed Brown 
4099371c9d4SSatish Balay static void THIViscosity(THI thi, PetscReal gam, PetscReal *eta, PetscReal *deta) {
410c4762a1bSJed Brown   PetscReal Bd2, eps, exponent;
411c4762a1bSJed Brown   if (thi->viscosity.Bd2 == 0) {
412c4762a1bSJed Brown     Units           units   = thi->units;
4139371c9d4SSatish Balay     const PetscReal n       = thi->viscosity.glen_n,                                  /* Glen exponent */
414c4762a1bSJed Brown       p                     = 1. + 1. / n,                                            /* for Stokes */
415c4762a1bSJed Brown       A                     = 1.e-16 * PetscPowReal(units->Pascal, -n) / units->year, /* softness parameter (Pa^{-n}/s) */
416c4762a1bSJed Brown       B                     = PetscPowReal(A, -1. / n);                               /* hardness parameter */
417c4762a1bSJed Brown     thi->viscosity.Bd2      = B / 2;
418c4762a1bSJed Brown     thi->viscosity.exponent = (p - 2) / 2;
419c4762a1bSJed Brown     thi->viscosity.eps      = 0.5 * PetscSqr(1e-5 / units->year);
420c4762a1bSJed Brown   }
421c4762a1bSJed Brown   Bd2      = thi->viscosity.Bd2;
422c4762a1bSJed Brown   exponent = thi->viscosity.exponent;
423c4762a1bSJed Brown   eps      = thi->viscosity.eps;
424c4762a1bSJed Brown   *eta     = Bd2 * PetscPowReal(eps + gam, exponent);
425c4762a1bSJed Brown   *deta    = exponent * (*eta) / (eps + gam);
426c4762a1bSJed Brown }
427c4762a1bSJed Brown 
4289371c9d4SSatish Balay static void THIErosion(THI thi, const Node *vel, PetscScalar *erate, Node *derate) {
4299371c9d4SSatish Balay   const PetscScalar magref2 = 1.e-10 + (PetscSqr(vel->u) + PetscSqr(vel->v)) / PetscSqr(thi->erosion.refvel), rate = -thi->erosion.rate * PetscPowScalar(magref2, 0.5 * thi->erosion.exponent);
430c4762a1bSJed Brown   if (erate) *erate = rate;
431c4762a1bSJed Brown   if (derate) {
432c4762a1bSJed Brown     if (thi->erosion.exponent == 1) {
433c4762a1bSJed Brown       derate->u = 0;
434c4762a1bSJed Brown       derate->v = 0;
435c4762a1bSJed Brown     } else {
436c4762a1bSJed Brown       derate->u = 0.5 * thi->erosion.exponent * rate / magref2 * 2. * vel->u / PetscSqr(thi->erosion.refvel);
437c4762a1bSJed Brown       derate->v = 0.5 * thi->erosion.exponent * rate / magref2 * 2. * vel->v / PetscSqr(thi->erosion.refvel);
438c4762a1bSJed Brown     }
439c4762a1bSJed Brown   }
440c4762a1bSJed Brown }
441c4762a1bSJed Brown 
4429371c9d4SSatish Balay static void RangeUpdate(PetscReal *min, PetscReal *max, PetscReal x) {
443c4762a1bSJed Brown   if (x < *min) *min = x;
444c4762a1bSJed Brown   if (x > *max) *max = x;
445c4762a1bSJed Brown }
446c4762a1bSJed Brown 
4479371c9d4SSatish Balay static void PRangeClear(PRange *p) {
448c4762a1bSJed Brown   p->cmin = p->min = 1e100;
449c4762a1bSJed Brown   p->cmax = p->max = -1e100;
450c4762a1bSJed Brown }
451c4762a1bSJed Brown 
4529371c9d4SSatish Balay static PetscErrorCode PRangeMinMax(PRange *p, PetscReal min, PetscReal max) {
453c4762a1bSJed Brown   PetscFunctionBeginUser;
454c4762a1bSJed Brown   p->cmin = min;
455c4762a1bSJed Brown   p->cmax = max;
456c4762a1bSJed Brown   if (min < p->min) p->min = min;
457c4762a1bSJed Brown   if (max > p->max) p->max = max;
458c4762a1bSJed Brown   PetscFunctionReturn(0);
459c4762a1bSJed Brown }
460c4762a1bSJed Brown 
4619371c9d4SSatish Balay static PetscErrorCode THIDestroy(THI *thi) {
462c4762a1bSJed Brown   PetscFunctionBeginUser;
463c4762a1bSJed Brown   if (--((PetscObject)(*thi))->refct > 0) PetscFunctionReturn(0);
4649566063dSJacob Faibussowitsch   PetscCall(PetscFree((*thi)->units));
4659566063dSJacob Faibussowitsch   PetscCall(PetscFree((*thi)->mattype));
4669566063dSJacob Faibussowitsch   PetscCall(PetscFree((*thi)->monitor_basename));
4679566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(thi));
468c4762a1bSJed Brown   PetscFunctionReturn(0);
469c4762a1bSJed Brown }
470c4762a1bSJed Brown 
4719371c9d4SSatish Balay static PetscErrorCode THICreate(MPI_Comm comm, THI *inthi) {
472c4762a1bSJed Brown   static PetscBool registered = PETSC_FALSE;
473c4762a1bSJed Brown   THI              thi;
474c4762a1bSJed Brown   Units            units;
475c4762a1bSJed Brown   char             monitor_basename[PETSC_MAX_PATH_LEN] = "thi-";
476c4762a1bSJed Brown   PetscErrorCode   ierr;
477c4762a1bSJed Brown 
478c4762a1bSJed Brown   PetscFunctionBeginUser;
479c4762a1bSJed Brown   *inthi = 0;
480c4762a1bSJed Brown   if (!registered) {
4819566063dSJacob Faibussowitsch     PetscCall(PetscClassIdRegister("Toy Hydrostatic Ice", &THI_CLASSID));
482c4762a1bSJed Brown     registered = PETSC_TRUE;
483c4762a1bSJed Brown   }
4849566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(thi, THI_CLASSID, "THI", "Toy Hydrostatic Ice", "THI", comm, THIDestroy, 0));
485c4762a1bSJed Brown 
4869566063dSJacob Faibussowitsch   PetscCall(PetscNew(&thi->units));
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   units           = thi->units;
489c4762a1bSJed Brown   units->meter    = 1e-2;
490c4762a1bSJed Brown   units->second   = 1e-7;
491c4762a1bSJed Brown   units->kilogram = 1e-12;
492c4762a1bSJed Brown 
493d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Scaled units options", "");
494c4762a1bSJed Brown   {
4959566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-units_meter", "1 meter in scaled length units", "", units->meter, &units->meter, NULL));
4969566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-units_second", "1 second in scaled time units", "", units->second, &units->second, NULL));
4979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-units_kilogram", "1 kilogram in scaled mass units", "", units->kilogram, &units->kilogram, NULL));
498c4762a1bSJed Brown   }
499d0609cedSBarry Smith   PetscOptionsEnd();
500c4762a1bSJed Brown   units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second));
501c4762a1bSJed Brown   units->year   = 31556926. * units->second, /* seconds per year */
502c4762a1bSJed Brown 
503c4762a1bSJed Brown     thi->Lx            = 10.e3;
504c4762a1bSJed Brown   thi->Ly              = 10.e3;
505c4762a1bSJed Brown   thi->Lz              = 1000;
506c4762a1bSJed Brown   thi->nlevels         = 1;
507c4762a1bSJed Brown   thi->dirichlet_scale = 1;
508c4762a1bSJed Brown   thi->verbose         = PETSC_FALSE;
509c4762a1bSJed Brown 
510c4762a1bSJed Brown   thi->viscosity.glen_n = 3.;
511c4762a1bSJed Brown   thi->erosion.rate     = 1e-3; /* m/a */
512c4762a1bSJed Brown   thi->erosion.exponent = 1.;
513c4762a1bSJed Brown   thi->erosion.refvel   = 1.; /* m/a */
514c4762a1bSJed Brown 
515d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Toy Hydrostatic Ice options", "");
516c4762a1bSJed Brown   {
517c4762a1bSJed Brown     QuadratureType quad       = QUAD_GAUSS;
518c4762a1bSJed Brown     char           homexp[]   = "A";
519c4762a1bSJed Brown     char           mtype[256] = MATSBAIJ;
520c4762a1bSJed Brown     PetscReal      L, m = 1.0;
521c4762a1bSJed Brown     PetscBool      flg;
522c4762a1bSJed Brown     L = thi->Lx;
5239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_L", "Domain size (m)", "", L, &L, &flg));
524c4762a1bSJed Brown     if (flg) thi->Lx = thi->Ly = L;
5259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_Lx", "X Domain size (m)", "", thi->Lx, &thi->Lx, NULL));
5269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_Ly", "Y Domain size (m)", "", thi->Ly, &thi->Ly, NULL));
5279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_Lz", "Z Domain size (m)", "", thi->Lz, &thi->Lz, NULL));
5289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-thi_hom", "ISMIP-HOM experiment (A or C)", "", homexp, homexp, sizeof(homexp), NULL));
529c4762a1bSJed Brown     switch (homexp[0] = toupper(homexp[0])) {
530c4762a1bSJed Brown     case 'A':
531c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_A;
532c4762a1bSJed Brown       thi->no_slip    = PETSC_TRUE;
533c4762a1bSJed Brown       thi->alpha      = 0.5;
534c4762a1bSJed Brown       break;
535c4762a1bSJed Brown     case 'C':
536c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_C;
537c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
538c4762a1bSJed Brown       thi->alpha      = 0.1;
539c4762a1bSJed Brown       break;
540c4762a1bSJed Brown     case 'F':
541c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_F;
542c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
543c4762a1bSJed Brown       thi->alpha      = 0.5;
544c4762a1bSJed Brown       break;
545c4762a1bSJed Brown     case 'X':
546c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_X;
547c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
548c4762a1bSJed Brown       thi->alpha      = 0.3;
549c4762a1bSJed Brown       break;
550c4762a1bSJed Brown     case 'Y':
551c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_Y;
552c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
553c4762a1bSJed Brown       thi->alpha      = 0.5;
554c4762a1bSJed Brown       break;
555c4762a1bSJed Brown     case 'Z':
556c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_Z;
557c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
558c4762a1bSJed Brown       thi->alpha      = 0.5;
559c4762a1bSJed Brown       break;
5609371c9d4SSatish Balay     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "HOM experiment '%c' not implemented", homexp[0]);
561c4762a1bSJed Brown     }
5629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-thi_quadrature", "Quadrature to use for 3D elements", "", QuadratureTypes, (PetscEnum)quad, (PetscEnum *)&quad, NULL));
563c4762a1bSJed Brown     switch (quad) {
564c4762a1bSJed Brown     case QUAD_GAUSS:
565c4762a1bSJed Brown       HexQInterp = HexQInterp_Gauss;
566c4762a1bSJed Brown       HexQDeriv  = HexQDeriv_Gauss;
567c4762a1bSJed Brown       break;
568c4762a1bSJed Brown     case QUAD_LOBATTO:
569c4762a1bSJed Brown       HexQInterp = HexQInterp_Lobatto;
570c4762a1bSJed Brown       HexQDeriv  = HexQDeriv_Lobatto;
571c4762a1bSJed Brown       break;
572c4762a1bSJed Brown     }
5739566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_alpha", "Bed angle (degrees)", "", thi->alpha, &thi->alpha, NULL));
5749566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_viscosity_glen_n", "Exponent in Glen flow law, 1=linear, infty=ideal plastic", NULL, thi->viscosity.glen_n, &thi->viscosity.glen_n, NULL));
5759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_friction_m", "Friction exponent, 0=Coulomb, 1=Navier", "", m, &m, NULL));
576c4762a1bSJed Brown     thi->friction.exponent = (m - 1) / 2;
5779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_erosion_rate", "Rate of erosion relative to sliding velocity at reference velocity (m/a)", NULL, thi->erosion.rate, &thi->erosion.rate, NULL));
5789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_erosion_exponent", "Power of sliding velocity appearing in erosion relation", NULL, thi->erosion.exponent, &thi->erosion.exponent, NULL));
5799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_erosion_refvel", "Reference sliding velocity for erosion (m/a)", NULL, thi->erosion.refvel, &thi->erosion.refvel, NULL));
580c4762a1bSJed Brown     thi->erosion.rate *= units->meter / units->year;
581c4762a1bSJed Brown     thi->erosion.refvel *= units->meter / units->year;
5829566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_dirichlet_scale", "Scale Dirichlet boundary conditions by this factor", "", thi->dirichlet_scale, &thi->dirichlet_scale, NULL));
5839566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_ssa_friction_scale", "Scale slip boundary conditions by this factor in SSA (2D) assembly", "", thi->ssa_friction_scale, &thi->ssa_friction_scale, NULL));
5849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-thi_inertia", "Coefficient of accelaration term in velocity system, physical is almost zero", NULL, thi->inertia, &thi->inertia, NULL));
5859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-thi_nlevels", "Number of levels of refinement", "", thi->nlevels, &thi->nlevels, NULL));
5869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-thi_mat_type", "Matrix type", "MatSetType", MatList, mtype, (char *)mtype, sizeof(mtype), NULL));
5879566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(mtype, &thi->mattype));
5889566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-thi_verbose", "Enable verbose output (like matrix sizes and statistics)", "", thi->verbose, &thi->verbose, NULL));
5899566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-thi_monitor", "Basename to write state files to", NULL, monitor_basename, monitor_basename, sizeof(monitor_basename), &flg));
590c4762a1bSJed Brown     if (flg) {
5919566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(monitor_basename, &thi->monitor_basename));
592c4762a1bSJed Brown       thi->monitor_interval = 1;
5939566063dSJacob Faibussowitsch       PetscCall(PetscOptionsInt("-thi_monitor_interval", "Frequency at which to write state files", NULL, thi->monitor_interval, &thi->monitor_interval, NULL));
594c4762a1bSJed Brown     }
595c4762a1bSJed Brown   }
596d0609cedSBarry Smith   PetscOptionsEnd();
597c4762a1bSJed Brown 
598c4762a1bSJed Brown   /* dimensionalize */
599c4762a1bSJed Brown   thi->Lx *= units->meter;
600c4762a1bSJed Brown   thi->Ly *= units->meter;
601c4762a1bSJed Brown   thi->Lz *= units->meter;
602c4762a1bSJed Brown   thi->alpha *= PETSC_PI / 180;
603c4762a1bSJed Brown 
604c4762a1bSJed Brown   PRangeClear(&thi->eta);
605c4762a1bSJed Brown   PRangeClear(&thi->beta2);
606c4762a1bSJed Brown 
607c4762a1bSJed Brown   {
6089371c9d4SSatish Balay     PetscReal u = 1000 * units->meter / (3e7 * units->second), gradu = u / (100 * units->meter), eta, deta, rho = 910 * units->kilogram / PetscPowRealInt(units->meter, 3), grav = 9.81 * units->meter / PetscSqr(units->second),
609c4762a1bSJed Brown               driving = rho * grav * PetscSinReal(thi->alpha) * 1000 * units->meter;
610c4762a1bSJed Brown     THIViscosity(thi, 0.5 * gradu * gradu, &eta, &deta);
611c4762a1bSJed Brown     thi->rhog = rho * grav;
612c4762a1bSJed Brown     if (thi->verbose) {
61363a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi), "Units: meter %8.2g  second %8.2g  kg %8.2g  Pa %8.2g\n", (double)units->meter, (double)units->second, (double)units->kilogram, (double)units->Pascal));
61463a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi), "Domain (%6.2g,%6.2g,%6.2g), pressure %8.2g, driving stress %8.2g\n", (double)thi->Lx, (double)thi->Ly, (double)thi->Lz, (double)(rho * grav * 1e3 * units->meter), (double)driving));
61563a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi), "Large velocity 1km/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n", (double)u, (double)gradu, (double)eta, (double)(2 * eta * gradu, 2 * eta * gradu / driving)));
616c4762a1bSJed Brown       THIViscosity(thi, 0.5 * PetscSqr(1e-3 * gradu), &eta, &deta);
61763a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi), "Small velocity 1m/a  %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n", (double)(1e-3 * u), (double)(1e-3 * gradu), (double)eta, (double)(2 * eta * 1e-3 * gradu, 2 * eta * 1e-3 * gradu / driving)));
618c4762a1bSJed Brown     }
619c4762a1bSJed Brown   }
620c4762a1bSJed Brown 
621c4762a1bSJed Brown   *inthi = thi;
622c4762a1bSJed Brown   PetscFunctionReturn(0);
623c4762a1bSJed Brown }
624c4762a1bSJed Brown 
625c4762a1bSJed Brown /* Our problem is periodic, but the domain has a mean slope of alpha so the bed does not line up between the upstream
626c4762a1bSJed Brown  * and downstream ends of the domain.  This function fixes the ghost values so that the domain appears truly periodic in
627c4762a1bSJed Brown  * the horizontal. */
6289371c9d4SSatish Balay static PetscErrorCode THIFixGhosts(THI thi, DM da3, DM da2, Vec X3, Vec X2) {
629c4762a1bSJed Brown   DMDALocalInfo info;
630c4762a1bSJed Brown   PrmNode     **x2;
631c4762a1bSJed Brown   PetscInt      i, j;
632c4762a1bSJed Brown 
633c4762a1bSJed Brown   PetscFunctionBeginUser;
6349566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da3, &info));
6359566063dSJacob Faibussowitsch   /* PetscCall(VecView(X2,PETSC_VIEWER_STDOUT_WORLD)); */
6369566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da2, X2, &x2));
637c4762a1bSJed Brown   for (i = info.gzs; i < info.gzs + info.gzm; i++) {
638c4762a1bSJed Brown     if (i > -1 && i < info.mz) continue;
6399371c9d4SSatish Balay     for (j = info.gys; j < info.gys + info.gym; j++) { x2[i][j].b += PetscSinReal(thi->alpha) * thi->Lx * (i < 0 ? 1.0 : -1.0); }
640c4762a1bSJed Brown   }
6419566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da2, X2, &x2));
6429566063dSJacob Faibussowitsch   /* PetscCall(VecView(X2,PETSC_VIEWER_STDOUT_WORLD)); */
643c4762a1bSJed Brown   PetscFunctionReturn(0);
644c4762a1bSJed Brown }
645c4762a1bSJed Brown 
6469371c9d4SSatish Balay static PetscErrorCode THIInitializePrm(THI thi, DM da2prm, PrmNode **p) {
647c4762a1bSJed Brown   PetscInt i, j, xs, xm, ys, ym, mx, my;
648c4762a1bSJed Brown 
649c4762a1bSJed Brown   PetscFunctionBeginUser;
6509566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da2prm, &ys, &xs, 0, &ym, &xm, 0));
6519566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da2prm, 0, &my, &mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
652c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
653c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
654c4762a1bSJed Brown       PetscReal xx = thi->Lx * i / mx, yy = thi->Ly * j / my;
655c4762a1bSJed Brown       thi->initialize(thi, xx, yy, &p[i][j]);
656c4762a1bSJed Brown     }
657c4762a1bSJed Brown   }
658c4762a1bSJed Brown   PetscFunctionReturn(0);
659c4762a1bSJed Brown }
660c4762a1bSJed Brown 
6619371c9d4SSatish Balay static PetscErrorCode THIInitial(THI thi, DM pack, Vec X) {
662c4762a1bSJed Brown   DM        da3, da2;
663c4762a1bSJed Brown   PetscInt  i, j, k, xs, xm, ys, ym, zs, zm, mx, my;
664c4762a1bSJed Brown   PetscReal hx, hy;
665c4762a1bSJed Brown   PrmNode **prm;
666c4762a1bSJed Brown   Node   ***x;
667c4762a1bSJed Brown   Vec       X3g, X2g, X2;
668c4762a1bSJed Brown 
669c4762a1bSJed Brown   PetscFunctionBeginUser;
6709566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntries(pack, &da3, &da2));
6719566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccess(pack, X, &X3g, &X2g));
6729566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da2, &X2));
673c4762a1bSJed Brown 
6749566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da3, 0, 0, &my, &mx, 0, 0, 0, 0, 0, 0, 0, 0, 0));
6759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da3, &zs, &ys, &xs, &zm, &ym, &xm));
6769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da3, X3g, &x));
6779566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da2, X2, &prm));
678c4762a1bSJed Brown 
6799566063dSJacob Faibussowitsch   PetscCall(THIInitializePrm(thi, da2, prm));
680c4762a1bSJed Brown 
681c4762a1bSJed Brown   hx = thi->Lx / mx;
682c4762a1bSJed Brown   hy = thi->Ly / my;
683c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
684c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
685c4762a1bSJed Brown       for (k = zs; k < zs + zm; k++) {
6869371c9d4SSatish Balay         const PetscScalar zm1 = zm - 1, drivingx = thi->rhog * (prm[i + 1][j].b + prm[i + 1][j].h - prm[i - 1][j].b - prm[i - 1][j].h) / (2 * hx), drivingy = thi->rhog * (prm[i][j + 1].b + prm[i][j + 1].h - prm[i][j - 1].b - prm[i][j - 1].h) / (2 * hy);
687c4762a1bSJed Brown         x[i][j][k].u = 0. * drivingx * prm[i][j].h * (PetscScalar)k / zm1;
688c4762a1bSJed Brown         x[i][j][k].v = 0. * drivingy * prm[i][j].h * (PetscScalar)k / zm1;
689c4762a1bSJed Brown       }
690c4762a1bSJed Brown     }
691c4762a1bSJed Brown   }
692c4762a1bSJed Brown 
6939566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da3, X3g, &x));
6949566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da2, X2, &prm));
695c4762a1bSJed Brown 
6969566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da2, X2, INSERT_VALUES, X2g));
6979566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da2, X2, INSERT_VALUES, X2g));
6989566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da2, &X2));
699c4762a1bSJed Brown 
7009566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccess(pack, X, &X3g, &X2g));
701c4762a1bSJed Brown   PetscFunctionReturn(0);
702c4762a1bSJed Brown }
703c4762a1bSJed Brown 
7049371c9d4SSatish Balay static void PointwiseNonlinearity(THI thi, const Node n[restrict 8], const PetscReal phi[restrict 3], PetscReal dphi[restrict 8][3], PetscScalar *restrict u, PetscScalar *restrict v, PetscScalar du[restrict 3], PetscScalar dv[restrict 3], PetscReal *eta, PetscReal *deta) {
705c4762a1bSJed Brown   PetscInt    l, ll;
706c4762a1bSJed Brown   PetscScalar gam;
707c4762a1bSJed Brown 
708c4762a1bSJed Brown   du[0] = du[1] = du[2] = 0;
709c4762a1bSJed Brown   dv[0] = dv[1] = dv[2] = 0;
710c4762a1bSJed Brown   *u                    = 0;
711c4762a1bSJed Brown   *v                    = 0;
712c4762a1bSJed Brown   for (l = 0; l < 8; l++) {
713c4762a1bSJed Brown     *u += phi[l] * n[l].u;
714c4762a1bSJed Brown     *v += phi[l] * n[l].v;
715c4762a1bSJed Brown     for (ll = 0; ll < 3; ll++) {
716c4762a1bSJed Brown       du[ll] += dphi[l][ll] * n[l].u;
717c4762a1bSJed Brown       dv[ll] += dphi[l][ll] * n[l].v;
718c4762a1bSJed Brown     }
719c4762a1bSJed Brown   }
720c4762a1bSJed Brown   gam = Sqr(du[0]) + Sqr(dv[1]) + du[0] * dv[1] + 0.25 * Sqr(du[1] + dv[0]) + 0.25 * Sqr(du[2]) + 0.25 * Sqr(dv[2]);
721c4762a1bSJed Brown   THIViscosity(thi, PetscRealPart(gam), eta, deta);
722c4762a1bSJed Brown }
723c4762a1bSJed Brown 
7249371c9d4SSatish Balay static PetscErrorCode THIFunctionLocal_3D(DMDALocalInfo *info, const Node ***x, const PrmNode **prm, const Node ***xdot, Node ***f, THI thi) {
725c4762a1bSJed Brown   PetscInt  xs, ys, xm, ym, zm, i, j, k, q, l;
726c4762a1bSJed Brown   PetscReal hx, hy, etamin, etamax, beta2min, beta2max;
727c4762a1bSJed Brown 
728c4762a1bSJed Brown   PetscFunctionBeginUser;
729c4762a1bSJed Brown   xs = info->zs;
730c4762a1bSJed Brown   ys = info->ys;
731c4762a1bSJed Brown   xm = info->zm;
732c4762a1bSJed Brown   ym = info->ym;
733c4762a1bSJed Brown   zm = info->xm;
734c4762a1bSJed Brown   hx = thi->Lx / info->mz;
735c4762a1bSJed Brown   hy = thi->Ly / info->my;
736c4762a1bSJed Brown 
737c4762a1bSJed Brown   etamin   = 1e100;
738c4762a1bSJed Brown   etamax   = 0;
739c4762a1bSJed Brown   beta2min = 1e100;
740c4762a1bSJed Brown   beta2max = 0;
741c4762a1bSJed Brown 
742c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
743c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
744c4762a1bSJed Brown       PrmNode pn[4], dpn[4][2];
745c4762a1bSJed Brown       QuadExtract(prm, i, j, pn);
7469566063dSJacob Faibussowitsch       PetscCall(QuadComputeGrad4(QuadQDeriv, hx, hy, pn, dpn));
747c4762a1bSJed Brown       for (k = 0; k < zm - 1; k++) {
748c4762a1bSJed Brown         PetscInt  ls = 0;
749c4762a1bSJed Brown         Node      n[8], ndot[8], *fn[8];
750c4762a1bSJed Brown         PetscReal zn[8], etabase = 0;
7512f613bf5SBarry Smith 
752c4762a1bSJed Brown         PrmHexGetZ(pn, k, zm, zn);
753c4762a1bSJed Brown         HexExtract(x, i, j, k, n);
7542f613bf5SBarry Smith         HexExtract(xdot, i, j, k, ndot);
755c4762a1bSJed Brown         HexExtractRef(f, i, j, k, fn);
756c4762a1bSJed Brown         if (thi->no_slip && k == 0) {
757c4762a1bSJed Brown           for (l = 0; l < 4; l++) n[l].u = n[l].v = 0;
758c4762a1bSJed Brown           /* The first 4 basis functions lie on the bottom layer, so their contribution is exactly 0, hence we can skip them */
759c4762a1bSJed Brown           ls = 4;
760c4762a1bSJed Brown         }
761c4762a1bSJed Brown         for (q = 0; q < 8; q++) {
762c4762a1bSJed Brown           PetscReal   dz[3], phi[8], dphi[8][3], jw, eta, deta;
763c4762a1bSJed Brown           PetscScalar du[3], dv[3], u, v, udot = 0, vdot = 0;
764c4762a1bSJed Brown           for (l = ls; l < 8; l++) {
765c4762a1bSJed Brown             udot += HexQInterp[q][l] * ndot[l].u;
766c4762a1bSJed Brown             vdot += HexQInterp[q][l] * ndot[l].v;
767c4762a1bSJed Brown           }
768c4762a1bSJed Brown           HexGrad(HexQDeriv[q], zn, dz);
769c4762a1bSJed Brown           HexComputeGeometry(q, hx, hy, dz, phi, dphi, &jw);
770c4762a1bSJed Brown           PointwiseNonlinearity(thi, n, phi, dphi, &u, &v, du, dv, &eta, &deta);
771c4762a1bSJed Brown           jw /= thi->rhog; /* scales residuals to be O(1) */
772c4762a1bSJed Brown           if (q == 0) etabase = eta;
773c4762a1bSJed Brown           RangeUpdate(&etamin, &etamax, eta);
774c4762a1bSJed Brown           for (l = ls; l < 8; l++) { /* test functions */
775c4762a1bSJed Brown             const PetscScalar ds[2] = {dpn[q % 4][0].h + dpn[q % 4][0].b, dpn[q % 4][1].h + dpn[q % 4][1].b};
776c4762a1bSJed Brown             const PetscReal   pp = phi[l], *dp = dphi[l];
777c4762a1bSJed Brown             fn[l]->u += dp[0] * jw * eta * (4. * du[0] + 2. * dv[1]) + dp[1] * jw * eta * (du[1] + dv[0]) + dp[2] * jw * eta * du[2] + pp * jw * thi->rhog * ds[0];
778c4762a1bSJed Brown             fn[l]->v += dp[1] * jw * eta * (2. * du[0] + 4. * dv[1]) + dp[0] * jw * eta * (du[1] + dv[0]) + dp[2] * jw * eta * dv[2] + pp * jw * thi->rhog * ds[1];
779c4762a1bSJed Brown             fn[l]->u += pp * jw * udot * thi->inertia * pp;
780c4762a1bSJed Brown             fn[l]->v += pp * jw * vdot * thi->inertia * pp;
781c4762a1bSJed Brown           }
782c4762a1bSJed Brown         }
783c4762a1bSJed Brown         if (k == 0) { /* we are on a bottom face */
784c4762a1bSJed Brown           if (thi->no_slip) {
785c4762a1bSJed Brown             /* Note: Non-Galerkin coarse grid operators are very sensitive to the scaling of Dirichlet boundary
786c4762a1bSJed Brown             * conditions.  After shenanigans above, etabase contains the effective viscosity at the closest quadrature
787c4762a1bSJed Brown             * point to the bed.  We want the diagonal entry in the Dirichlet condition to have similar magnitude to the
788c4762a1bSJed Brown             * diagonal entry corresponding to the adjacent node.  The fundamental scaling of the viscous part is in
789c4762a1bSJed Brown             * diagu, diagv below.  This scaling is easy to recognize by considering the finite difference operator after
790c4762a1bSJed Brown             * scaling by element size.  The no-slip Dirichlet condition is scaled by this factor, and also in the
791c4762a1bSJed Brown             * assembled matrix (see the similar block in THIJacobianLocal).
792c4762a1bSJed Brown             *
793c4762a1bSJed Brown             * Note that the residual at this Dirichlet node is linear in the state at this node, but also depends
794c4762a1bSJed Brown             * (nonlinearly in general) on the neighboring interior nodes through the local viscosity.  This will make
795c4762a1bSJed Brown             * a matrix-free Jacobian have extra entries in the corresponding row.  We assemble only the diagonal part,
796c4762a1bSJed Brown             * so the solution will exactly satisfy the boundary condition after the first linear iteration.
797c4762a1bSJed Brown             */
798c4762a1bSJed Brown             const PetscReal   hz    = PetscRealPart(pn[0].h) / (zm - 1.);
799c4762a1bSJed Brown             const PetscScalar diagu = 2 * etabase / thi->rhog * (hx * hy / hz + hx * hz / hy + 4 * hy * hz / hx), diagv = 2 * etabase / thi->rhog * (hx * hy / hz + 4 * hx * hz / hy + hy * hz / hx);
800c4762a1bSJed Brown             fn[0]->u = thi->dirichlet_scale * diagu * x[i][j][k].u;
801c4762a1bSJed Brown             fn[0]->v = thi->dirichlet_scale * diagv * x[i][j][k].v;
802c4762a1bSJed Brown           } else {                    /* Integrate over bottom face to apply boundary condition */
803c4762a1bSJed Brown             for (q = 0; q < 4; q++) { /* We remove the explicit scaling of the residual by 1/rhog because beta2 already has that scaling to be O(1) */
804c4762a1bSJed Brown               const PetscReal jw = 0.25 * hx * hy, *phi = QuadQInterp[q];
805c4762a1bSJed Brown               PetscScalar     u = 0, v = 0, rbeta2 = 0;
806c4762a1bSJed Brown               PetscReal       beta2, dbeta2;
807c4762a1bSJed Brown               for (l = 0; l < 4; l++) {
808c4762a1bSJed Brown                 u += phi[l] * n[l].u;
809c4762a1bSJed Brown                 v += phi[l] * n[l].v;
810c4762a1bSJed Brown                 rbeta2 += phi[l] * pn[l].beta2;
811c4762a1bSJed Brown               }
812c4762a1bSJed Brown               THIFriction(thi, PetscRealPart(rbeta2), PetscRealPart(u * u + v * v) / 2, &beta2, &dbeta2);
813c4762a1bSJed Brown               RangeUpdate(&beta2min, &beta2max, beta2);
814c4762a1bSJed Brown               for (l = 0; l < 4; l++) {
815c4762a1bSJed Brown                 const PetscReal pp = phi[l];
816c4762a1bSJed Brown                 fn[ls + l]->u += pp * jw * beta2 * u;
817c4762a1bSJed Brown                 fn[ls + l]->v += pp * jw * beta2 * v;
818c4762a1bSJed Brown               }
819c4762a1bSJed Brown             }
820c4762a1bSJed Brown           }
821c4762a1bSJed Brown         }
822c4762a1bSJed Brown       }
823c4762a1bSJed Brown     }
824c4762a1bSJed Brown   }
825c4762a1bSJed Brown 
8269566063dSJacob Faibussowitsch   PetscCall(PRangeMinMax(&thi->eta, etamin, etamax));
8279566063dSJacob Faibussowitsch   PetscCall(PRangeMinMax(&thi->beta2, beta2min, beta2max));
828c4762a1bSJed Brown   PetscFunctionReturn(0);
829c4762a1bSJed Brown }
830c4762a1bSJed Brown 
8319371c9d4SSatish Balay static PetscErrorCode THIFunctionLocal_2D(DMDALocalInfo *info, const Node ***x, const PrmNode **prm, const PrmNode **prmdot, PrmNode **f, THI thi) {
832c4762a1bSJed Brown   PetscInt xs, ys, xm, ym, zm, i, j, k;
833c4762a1bSJed Brown 
834c4762a1bSJed Brown   PetscFunctionBeginUser;
835c4762a1bSJed Brown   xs = info->zs;
836c4762a1bSJed Brown   ys = info->ys;
837c4762a1bSJed Brown   xm = info->zm;
838c4762a1bSJed Brown   ym = info->ym;
839c4762a1bSJed Brown   zm = info->xm;
840c4762a1bSJed Brown 
841c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
842c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
843c4762a1bSJed Brown       PetscScalar div = 0, erate, h[8];
844c4762a1bSJed Brown       PrmNodeGetFaceMeasure(prm, i, j, h);
845c4762a1bSJed Brown       for (k = 0; k < zm; k++) {
846c4762a1bSJed Brown         PetscScalar weight = (k == 0 || k == zm - 1) ? 0.5 / (zm - 1) : 1.0 / (zm - 1);
847c4762a1bSJed Brown         if (0) { /* centered flux */
8489371c9d4SSatish Balay           div += (-weight * h[0] * StaggeredMidpoint2D(x[i][j][k].u, x[i - 1][j][k].u, x[i - 1][j - 1][k].u, x[i][j - 1][k].u) - weight * h[1] * StaggeredMidpoint2D(x[i][j][k].u, x[i - 1][j][k].u, x[i - 1][j + 1][k].u, x[i][j + 1][k].u) +
8499371c9d4SSatish Balay                   weight * h[2] * StaggeredMidpoint2D(x[i][j][k].u, x[i + 1][j][k].u, x[i + 1][j + 1][k].u, x[i][j + 1][k].u) + weight * h[3] * StaggeredMidpoint2D(x[i][j][k].u, x[i + 1][j][k].u, x[i + 1][j - 1][k].u, x[i][j - 1][k].u) -
8509371c9d4SSatish Balay                   weight * h[4] * StaggeredMidpoint2D(x[i][j][k].v, x[i][j - 1][k].v, x[i + 1][j - 1][k].v, x[i + 1][j][k].v) - weight * h[5] * StaggeredMidpoint2D(x[i][j][k].v, x[i][j - 1][k].v, x[i - 1][j - 1][k].v, x[i - 1][j][k].v) +
8519371c9d4SSatish Balay                   weight * h[6] * StaggeredMidpoint2D(x[i][j][k].v, x[i][j + 1][k].v, x[i - 1][j + 1][k].v, x[i - 1][j][k].v) + weight * h[7] * StaggeredMidpoint2D(x[i][j][k].v, x[i][j + 1][k].v, x[i + 1][j + 1][k].v, x[i + 1][j][k].v));
852c4762a1bSJed Brown         } else { /* Upwind flux */
8539371c9d4SSatish Balay           div += weight * (-UpwindFluxXW(x, prm, h, i, j, k, 1) - UpwindFluxXW(x, prm, h, i, j, k, -1) + UpwindFluxXE(x, prm, h, i, j, k, 1) + UpwindFluxXE(x, prm, h, i, j, k, -1) - UpwindFluxYS(x, prm, h, i, j, k, 1) - UpwindFluxYS(x, prm, h, i, j, k, -1) + UpwindFluxYN(x, prm, h, i, j, k, 1) + UpwindFluxYN(x, prm, h, i, j, k, -1));
854c4762a1bSJed Brown         }
855c4762a1bSJed Brown       }
856c4762a1bSJed Brown       /* printf("div[%d][%d] %g\n",i,j,div); */
857c4762a1bSJed Brown       THIErosion(thi, &x[i][j][0], &erate, NULL);
858c4762a1bSJed Brown       f[i][j].b     = prmdot[i][j].b - erate;
859c4762a1bSJed Brown       f[i][j].h     = prmdot[i][j].h + div;
860c4762a1bSJed Brown       f[i][j].beta2 = prmdot[i][j].beta2;
861c4762a1bSJed Brown     }
862c4762a1bSJed Brown   }
863c4762a1bSJed Brown   PetscFunctionReturn(0);
864c4762a1bSJed Brown }
865c4762a1bSJed Brown 
8669371c9d4SSatish Balay static PetscErrorCode THIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
867c4762a1bSJed Brown   THI             thi = (THI)ctx;
868c4762a1bSJed Brown   DM              pack, da3, da2;
869c4762a1bSJed Brown   Vec             X3, X2, Xdot3, Xdot2, F3, F2, F3g, F2g;
870c4762a1bSJed Brown   const Node   ***x3, ***xdot3;
871c4762a1bSJed Brown   const PrmNode **x2, **xdot2;
872c4762a1bSJed Brown   Node         ***f3;
873c4762a1bSJed Brown   PrmNode       **f2;
874c4762a1bSJed Brown   DMDALocalInfo   info3;
875c4762a1bSJed Brown 
876c4762a1bSJed Brown   PetscFunctionBeginUser;
8779566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &pack));
8789566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntries(pack, &da3, &da2));
8799566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da3, &info3));
8809566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(pack, &X3, &X2));
8819566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(pack, &Xdot3, &Xdot2));
8829566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(pack, X, X3, X2));
8839566063dSJacob Faibussowitsch   PetscCall(THIFixGhosts(thi, da3, da2, X3, X2));
8849566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(pack, Xdot, Xdot3, Xdot2));
885c4762a1bSJed Brown 
8869566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da3, &F3));
8879566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da2, &F2));
8889566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F3));
889c4762a1bSJed Brown 
8909566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da3, X3, &x3));
8919566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da2, X2, &x2));
8929566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da3, Xdot3, &xdot3));
8939566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da2, Xdot2, &xdot2));
8949566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da3, F3, &f3));
8959566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da2, F2, &f2));
896c4762a1bSJed Brown 
8979566063dSJacob Faibussowitsch   PetscCall(THIFunctionLocal_3D(&info3, x3, x2, xdot3, f3, thi));
8989566063dSJacob Faibussowitsch   PetscCall(THIFunctionLocal_2D(&info3, x3, x2, xdot2, f2, thi));
899c4762a1bSJed Brown 
9009566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da3, X3, &x3));
9019566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da2, X2, &x2));
9029566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da3, Xdot3, &xdot3));
9039566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da2, Xdot2, &xdot2));
9049566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da3, F3, &f3));
9059566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da2, F2, &f2));
906c4762a1bSJed Brown 
9079566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(pack, &X3, &X2));
9089566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(pack, &Xdot3, &Xdot2));
909c4762a1bSJed Brown 
9109566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
9119566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccess(pack, F, &F3g, &F2g));
9129566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da3, F3, ADD_VALUES, F3g));
9139566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da3, F3, ADD_VALUES, F3g));
9149566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da2, F2, INSERT_VALUES, F2g));
9159566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da2, F2, INSERT_VALUES, F2g));
916c4762a1bSJed Brown 
917c4762a1bSJed Brown   if (thi->verbose) {
918c4762a1bSJed Brown     PetscViewer viewer;
9199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)thi), &viewer));
9209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "3D_Velocity residual (bs=2):\n"));
9219566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
9229566063dSJacob Faibussowitsch     PetscCall(VecView(F3, viewer));
9239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
9249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "2D_Fields residual (bs=3):\n"));
9259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
9269566063dSJacob Faibussowitsch     PetscCall(VecView(F2, viewer));
9279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
928c4762a1bSJed Brown   }
929c4762a1bSJed Brown 
9309566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccess(pack, F, &F3g, &F2g));
931c4762a1bSJed Brown 
9329566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da3, &F3));
9339566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da2, &F2));
934c4762a1bSJed Brown   PetscFunctionReturn(0);
935c4762a1bSJed Brown }
936c4762a1bSJed Brown 
9379371c9d4SSatish Balay static PetscErrorCode THIMatrixStatistics(THI thi, Mat B, PetscViewer viewer) {
938c4762a1bSJed Brown   PetscReal   nrm;
939c4762a1bSJed Brown   PetscInt    m;
940c4762a1bSJed Brown   PetscMPIInt rank;
941c4762a1bSJed Brown 
942c4762a1bSJed Brown   PetscFunctionBeginUser;
9439566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_FROBENIUS, &nrm));
9449566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &m, 0));
9459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B), &rank));
946dd400576SPatrick Sanan   if (rank == 0) {
947c4762a1bSJed Brown     PetscScalar val0, val2;
9489566063dSJacob Faibussowitsch     PetscCall(MatGetValue(B, 0, 0, &val0));
9499566063dSJacob Faibussowitsch     PetscCall(MatGetValue(B, 2, 2, &val2));
9509371c9d4SSatish Balay     PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix dim %8" PetscInt_FMT "  norm %8.2e, (0,0) %8.2e  (2,2) %8.2e, eta [%8.2e,%8.2e] beta2 [%8.2e,%8.2e]\n", m, (double)nrm, (double)PetscRealPart(val0), (double)PetscRealPart(val2), (double)thi->eta.cmin,
9519371c9d4SSatish Balay                                      (double)thi->eta.cmax, (double)thi->beta2.cmin, (double)thi->beta2.cmax));
952c4762a1bSJed Brown   }
953c4762a1bSJed Brown   PetscFunctionReturn(0);
954c4762a1bSJed Brown }
955c4762a1bSJed Brown 
9569371c9d4SSatish Balay static PetscErrorCode THISurfaceStatistics(DM pack, Vec X, PetscReal *min, PetscReal *max, PetscReal *mean) {
957c4762a1bSJed Brown   DM          da3, da2;
958c4762a1bSJed Brown   Vec         X3, X2;
959c4762a1bSJed Brown   Node     ***x;
960c4762a1bSJed Brown   PetscInt    i, j, xs, ys, zs, xm, ym, zm, mx, my, mz;
961c4762a1bSJed Brown   PetscReal   umin = 1e100, umax = -1e100;
962c4762a1bSJed Brown   PetscScalar usum = 0.0, gusum;
963c4762a1bSJed Brown 
964c4762a1bSJed Brown   PetscFunctionBeginUser;
9659566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntries(pack, &da3, &da2));
9669566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccess(pack, X, &X3, &X2));
967c4762a1bSJed Brown   *min = *max = *mean = 0;
9689566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da3, 0, &mz, &my, &mx, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9699566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da3, &zs, &ys, &xs, &zm, &ym, &xm));
9703c633725SBarry Smith   PetscCheck(zs == 0 && zm == mz, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Unexpected decomposition");
9719566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da3, X3, &x));
972c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
973c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
974c4762a1bSJed Brown       PetscReal u = PetscRealPart(x[i][j][zm - 1].u);
975c4762a1bSJed Brown       RangeUpdate(&umin, &umax, u);
976c4762a1bSJed Brown       usum += u;
977c4762a1bSJed Brown     }
978c4762a1bSJed Brown   }
9799566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da3, X3, &x));
9809566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccess(pack, X, &X3, &X2));
981c4762a1bSJed Brown 
9829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&umin, min, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)da3)));
9839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&umax, max, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)da3)));
9849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&usum, &gusum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)da3)));
985c4762a1bSJed Brown   *mean = PetscRealPart(gusum) / (mx * my);
986c4762a1bSJed Brown   PetscFunctionReturn(0);
987c4762a1bSJed Brown }
988c4762a1bSJed Brown 
9899371c9d4SSatish Balay static PetscErrorCode THISolveStatistics(THI thi, TS ts, PetscInt coarsened, const char name[]) {
990c4762a1bSJed Brown   MPI_Comm comm;
991c4762a1bSJed Brown   DM       pack;
992c4762a1bSJed Brown   Vec      X, X3, X2;
993c4762a1bSJed Brown 
994c4762a1bSJed Brown   PetscFunctionBeginUser;
9959566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)thi, &comm));
9969566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &pack));
9979566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &X));
9989566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccess(pack, X, &X3, &X2));
9999566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Solution statistics after solve: %s\n", name));
1000c4762a1bSJed Brown   {
1001c4762a1bSJed Brown     PetscInt            its, lits;
1002c4762a1bSJed Brown     SNESConvergedReason reason;
1003c4762a1bSJed Brown     SNES                snes;
10049566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
10059566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &its));
10069566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes, &reason));
10079566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(snes, &lits));
100863a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "%s: Number of SNES iterations = %" PetscInt_FMT ", total linear iterations = %" PetscInt_FMT "\n", SNESConvergedReasons[reason], its, lits));
1009c4762a1bSJed Brown   }
1010c4762a1bSJed Brown   {
1011c4762a1bSJed Brown     PetscReal    nrm2, tmin[3] = {1e100, 1e100, 1e100}, tmax[3] = {-1e100, -1e100, -1e100}, min[3], max[3];
1012c4762a1bSJed Brown     PetscInt     i, j, m;
1013c4762a1bSJed Brown     PetscScalar *x;
10149566063dSJacob Faibussowitsch     PetscCall(VecNorm(X3, NORM_2, &nrm2));
10159566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(X3, &m));
10169566063dSJacob Faibussowitsch     PetscCall(VecGetArray(X3, &x));
1017c4762a1bSJed Brown     for (i = 0; i < m; i += 2) {
1018c4762a1bSJed Brown       PetscReal u = PetscRealPart(x[i]), v = PetscRealPart(x[i + 1]), c = PetscSqrtReal(u * u + v * v);
1019c4762a1bSJed Brown       tmin[0] = PetscMin(u, tmin[0]);
1020c4762a1bSJed Brown       tmin[1] = PetscMin(v, tmin[1]);
1021c4762a1bSJed Brown       tmin[2] = PetscMin(c, tmin[2]);
1022c4762a1bSJed Brown       tmax[0] = PetscMax(u, tmax[0]);
1023c4762a1bSJed Brown       tmax[1] = PetscMax(v, tmax[1]);
1024c4762a1bSJed Brown       tmax[2] = PetscMax(c, tmax[2]);
1025c4762a1bSJed Brown     }
10269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(X, &x));
10279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(tmin, min, 3, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)thi)));
10289566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(tmax, max, 3, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)thi)));
1029c4762a1bSJed Brown     /* Dimensionalize to meters/year */
1030c4762a1bSJed Brown     nrm2 *= thi->units->year / thi->units->meter;
1031c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1032c4762a1bSJed Brown       min[j] *= thi->units->year / thi->units->meter;
1033c4762a1bSJed Brown       max[j] *= thi->units->year / thi->units->meter;
1034c4762a1bSJed Brown     }
103563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "|X|_2 %g   u in [%g, %g]   v in [%g, %g]   c in [%g, %g] \n", (double)nrm2, (double)min[0], (double)max[0], (double)min[1], (double)max[1], (double)min[2], (double)max[2]));
1036c4762a1bSJed Brown     {
1037c4762a1bSJed Brown       PetscReal umin, umax, umean;
10389566063dSJacob Faibussowitsch       PetscCall(THISurfaceStatistics(pack, X, &umin, &umax, &umean));
1039c4762a1bSJed Brown       umin *= thi->units->year / thi->units->meter;
1040c4762a1bSJed Brown       umax *= thi->units->year / thi->units->meter;
1041c4762a1bSJed Brown       umean *= thi->units->year / thi->units->meter;
104263a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Surface statistics: u in [%12.6e, %12.6e] mean %12.6e\n", (double)umin, (double)umax, (double)umean));
1043c4762a1bSJed Brown     }
1044c4762a1bSJed Brown     /* These values stay nondimensional */
104563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Global eta range   [%g, %g], converged range [%g, %g]\n", (double)thi->eta.min, (double)thi->eta.max, (double)thi->eta.cmin, (double)thi->eta.cmax));
104663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Global beta2 range [%g, %g], converged range [%g, %g]\n", (double)thi->beta2.min, (double)thi->beta2.max, (double)thi->beta2.cmin, (double)thi->beta2.cmax));
1047c4762a1bSJed Brown   }
10489566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
10499566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccess(pack, X, &X3, &X2));
1050c4762a1bSJed Brown   PetscFunctionReturn(0);
1051c4762a1bSJed Brown }
1052c4762a1bSJed Brown 
10539371c9d4SSatish Balay static inline PetscInt DMDALocalIndex3D(DMDALocalInfo *info, PetscInt i, PetscInt j, PetscInt k) {
10549371c9d4SSatish Balay   return ((i - info->gzs) * info->gym + (j - info->gys)) * info->gxm + (k - info->gxs);
10559371c9d4SSatish Balay }
10569371c9d4SSatish Balay static inline PetscInt DMDALocalIndex2D(DMDALocalInfo *info, PetscInt i, PetscInt j) {
10579371c9d4SSatish Balay   return (i - info->gzs) * info->gym + (j - info->gys);
10589371c9d4SSatish Balay }
1059c4762a1bSJed Brown 
10609371c9d4SSatish Balay static PetscErrorCode THIJacobianLocal_Momentum(DMDALocalInfo *info, const Node ***x, const PrmNode **prm, Mat B, Mat Bcpl, THI thi) {
1061c4762a1bSJed Brown   PetscInt  xs, ys, xm, ym, zm, i, j, k, q, l, ll;
1062c4762a1bSJed Brown   PetscReal hx, hy;
1063c4762a1bSJed Brown 
1064c4762a1bSJed Brown   PetscFunctionBeginUser;
1065c4762a1bSJed Brown   xs = info->zs;
1066c4762a1bSJed Brown   ys = info->ys;
1067c4762a1bSJed Brown   xm = info->zm;
1068c4762a1bSJed Brown   ym = info->ym;
1069c4762a1bSJed Brown   zm = info->xm;
1070c4762a1bSJed Brown   hx = thi->Lx / info->mz;
1071c4762a1bSJed Brown   hy = thi->Ly / info->my;
1072c4762a1bSJed Brown 
1073c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
1074c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
1075c4762a1bSJed Brown       PrmNode pn[4], dpn[4][2];
1076c4762a1bSJed Brown       QuadExtract(prm, i, j, pn);
10779566063dSJacob Faibussowitsch       PetscCall(QuadComputeGrad4(QuadQDeriv, hx, hy, pn, dpn));
1078c4762a1bSJed Brown       for (k = 0; k < zm - 1; k++) {
1079c4762a1bSJed Brown         Node        n[8];
1080c4762a1bSJed Brown         PetscReal   zn[8], etabase = 0;
1081c4762a1bSJed Brown         PetscScalar Ke[8 * NODE_SIZE][8 * NODE_SIZE], Kcpl[8 * NODE_SIZE][4 * PRMNODE_SIZE];
1082c4762a1bSJed Brown         PetscInt    ls = 0;
1083c4762a1bSJed Brown 
1084c4762a1bSJed Brown         PrmHexGetZ(pn, k, zm, zn);
1085c4762a1bSJed Brown         HexExtract(x, i, j, k, n);
10869566063dSJacob Faibussowitsch         PetscCall(PetscMemzero(Ke, sizeof(Ke)));
10879566063dSJacob Faibussowitsch         PetscCall(PetscMemzero(Kcpl, sizeof(Kcpl)));
1088c4762a1bSJed Brown         if (thi->no_slip && k == 0) {
1089c4762a1bSJed Brown           for (l = 0; l < 4; l++) n[l].u = n[l].v = 0;
1090c4762a1bSJed Brown           ls = 4;
1091c4762a1bSJed Brown         }
1092c4762a1bSJed Brown         for (q = 0; q < 8; q++) {
1093c4762a1bSJed Brown           PetscReal   dz[3], phi[8], dphi[8][3], jw, eta, deta;
1094c4762a1bSJed Brown           PetscScalar du[3], dv[3], u, v;
1095c4762a1bSJed Brown           HexGrad(HexQDeriv[q], zn, dz);
1096c4762a1bSJed Brown           HexComputeGeometry(q, hx, hy, dz, phi, dphi, &jw);
1097c4762a1bSJed Brown           PointwiseNonlinearity(thi, n, phi, dphi, &u, &v, du, dv, &eta, &deta);
1098c4762a1bSJed Brown           jw /= thi->rhog; /* residuals are scaled by this factor */
1099c4762a1bSJed Brown           if (q == 0) etabase = eta;
1100c4762a1bSJed Brown           for (l = ls; l < 8; l++) { /* test functions */
1101c4762a1bSJed Brown             const PetscReal pp = phi[l], *restrict dp = dphi[l];
1102c4762a1bSJed Brown             for (ll = ls; ll < 8; ll++) { /* trial functions */
1103c4762a1bSJed Brown               const PetscReal *restrict dpl = dphi[ll];
1104c4762a1bSJed Brown               PetscScalar dgdu, dgdv;
1105c4762a1bSJed Brown               dgdu = 2. * du[0] * dpl[0] + dv[1] * dpl[0] + 0.5 * (du[1] + dv[0]) * dpl[1] + 0.5 * du[2] * dpl[2];
1106c4762a1bSJed Brown               dgdv = 2. * dv[1] * dpl[1] + du[0] * dpl[1] + 0.5 * (du[1] + dv[0]) * dpl[0] + 0.5 * dv[2] * dpl[2];
1107c4762a1bSJed Brown               /* Picard part */
1108c4762a1bSJed Brown               Ke[l * 2 + 0][ll * 2 + 0] += dp[0] * jw * eta * 4. * dpl[0] + dp[1] * jw * eta * dpl[1] + dp[2] * jw * eta * dpl[2];
1109c4762a1bSJed Brown               Ke[l * 2 + 0][ll * 2 + 1] += dp[0] * jw * eta * 2. * dpl[1] + dp[1] * jw * eta * dpl[0];
1110c4762a1bSJed Brown               Ke[l * 2 + 1][ll * 2 + 0] += dp[1] * jw * eta * 2. * dpl[0] + dp[0] * jw * eta * dpl[1];
1111c4762a1bSJed Brown               Ke[l * 2 + 1][ll * 2 + 1] += dp[1] * jw * eta * 4. * dpl[1] + dp[0] * jw * eta * dpl[0] + dp[2] * jw * eta * dpl[2];
1112c4762a1bSJed Brown               /* extra Newton terms */
1113c4762a1bSJed Brown               Ke[l * 2 + 0][ll * 2 + 0] += dp[0] * jw * deta * dgdu * (4. * du[0] + 2. * dv[1]) + dp[1] * jw * deta * dgdu * (du[1] + dv[0]) + dp[2] * jw * deta * dgdu * du[2];
1114c4762a1bSJed Brown               Ke[l * 2 + 0][ll * 2 + 1] += dp[0] * jw * deta * dgdv * (4. * du[0] + 2. * dv[1]) + dp[1] * jw * deta * dgdv * (du[1] + dv[0]) + dp[2] * jw * deta * dgdv * du[2];
1115c4762a1bSJed Brown               Ke[l * 2 + 1][ll * 2 + 0] += dp[1] * jw * deta * dgdu * (4. * dv[1] + 2. * du[0]) + dp[0] * jw * deta * dgdu * (du[1] + dv[0]) + dp[2] * jw * deta * dgdu * dv[2];
1116c4762a1bSJed Brown               Ke[l * 2 + 1][ll * 2 + 1] += dp[1] * jw * deta * dgdv * (4. * dv[1] + 2. * du[0]) + dp[0] * jw * deta * dgdv * (du[1] + dv[0]) + dp[2] * jw * deta * dgdv * dv[2];
1117c4762a1bSJed Brown               /* inertial part */
1118c4762a1bSJed Brown               Ke[l * 2 + 0][ll * 2 + 0] += pp * jw * thi->inertia * pp;
1119c4762a1bSJed Brown               Ke[l * 2 + 1][ll * 2 + 1] += pp * jw * thi->inertia * pp;
1120c4762a1bSJed Brown             }
1121c4762a1bSJed Brown             for (ll = 0; ll < 4; ll++) {                                                              /* Trial functions for surface/bed */
1122c4762a1bSJed Brown               const PetscReal dpl[] = {QuadQDeriv[q % 4][ll][0] / hx, QuadQDeriv[q % 4][ll][1] / hy}; /* surface = h + b */
1123c4762a1bSJed Brown               Kcpl[FieldIndex(Node, l, u)][FieldIndex(PrmNode, ll, h)] += pp * jw * thi->rhog * dpl[0];
1124c4762a1bSJed Brown               Kcpl[FieldIndex(Node, l, u)][FieldIndex(PrmNode, ll, b)] += pp * jw * thi->rhog * dpl[0];
1125c4762a1bSJed Brown               Kcpl[FieldIndex(Node, l, v)][FieldIndex(PrmNode, ll, h)] += pp * jw * thi->rhog * dpl[1];
1126c4762a1bSJed Brown               Kcpl[FieldIndex(Node, l, v)][FieldIndex(PrmNode, ll, b)] += pp * jw * thi->rhog * dpl[1];
1127c4762a1bSJed Brown             }
1128c4762a1bSJed Brown           }
1129c4762a1bSJed Brown         }
1130c4762a1bSJed Brown         if (k == 0) { /* on a bottom face */
1131c4762a1bSJed Brown           if (thi->no_slip) {
1132c4762a1bSJed Brown             const PetscReal   hz    = PetscRealPart(pn[0].h) / (zm - 1);
1133c4762a1bSJed Brown             const PetscScalar diagu = 2 * etabase / thi->rhog * (hx * hy / hz + hx * hz / hy + 4 * hy * hz / hx), diagv = 2 * etabase / thi->rhog * (hx * hy / hz + 4 * hx * hz / hy + hy * hz / hx);
1134c4762a1bSJed Brown             Ke[0][0] = thi->dirichlet_scale * diagu;
1135c4762a1bSJed Brown             Ke[0][1] = 0;
1136c4762a1bSJed Brown             Ke[1][0] = 0;
1137c4762a1bSJed Brown             Ke[1][1] = thi->dirichlet_scale * diagv;
1138c4762a1bSJed Brown           } else {
1139c4762a1bSJed Brown             for (q = 0; q < 4; q++) { /* We remove the explicit scaling by 1/rhog because beta2 already has that scaling to be O(1) */
1140c4762a1bSJed Brown               const PetscReal jw = 0.25 * hx * hy, *phi = QuadQInterp[q];
1141c4762a1bSJed Brown               PetscScalar     u = 0, v = 0, rbeta2 = 0;
1142c4762a1bSJed Brown               PetscReal       beta2, dbeta2;
1143c4762a1bSJed Brown               for (l = 0; l < 4; l++) {
1144c4762a1bSJed Brown                 u += phi[l] * n[l].u;
1145c4762a1bSJed Brown                 v += phi[l] * n[l].v;
1146c4762a1bSJed Brown                 rbeta2 += phi[l] * pn[l].beta2;
1147c4762a1bSJed Brown               }
1148c4762a1bSJed Brown               THIFriction(thi, PetscRealPart(rbeta2), PetscRealPart(u * u + v * v) / 2, &beta2, &dbeta2);
1149c4762a1bSJed Brown               for (l = 0; l < 4; l++) {
1150c4762a1bSJed Brown                 const PetscReal pp = phi[l];
1151c4762a1bSJed Brown                 for (ll = 0; ll < 4; ll++) {
1152c4762a1bSJed Brown                   const PetscReal ppl = phi[ll];
1153c4762a1bSJed Brown                   Ke[l * 2 + 0][ll * 2 + 0] += pp * jw * beta2 * ppl + pp * jw * dbeta2 * u * u * ppl;
1154c4762a1bSJed Brown                   Ke[l * 2 + 0][ll * 2 + 1] += pp * jw * dbeta2 * u * v * ppl;
1155c4762a1bSJed Brown                   Ke[l * 2 + 1][ll * 2 + 0] += pp * jw * dbeta2 * v * u * ppl;
1156c4762a1bSJed Brown                   Ke[l * 2 + 1][ll * 2 + 1] += pp * jw * beta2 * ppl + pp * jw * dbeta2 * v * v * ppl;
1157c4762a1bSJed Brown                 }
1158c4762a1bSJed Brown               }
1159c4762a1bSJed Brown             }
1160c4762a1bSJed Brown           }
1161c4762a1bSJed Brown         }
1162c4762a1bSJed Brown         {
11639371c9d4SSatish Balay           const PetscInt rc3blocked[8]                 = {DMDALocalIndex3D(info, i + 0, j + 0, k + 0), DMDALocalIndex3D(info, i + 1, j + 0, k + 0), DMDALocalIndex3D(info, i + 1, j + 1, k + 0), DMDALocalIndex3D(info, i + 0, j + 1, k + 0),
11649371c9d4SSatish Balay                                                           DMDALocalIndex3D(info, i + 0, j + 0, k + 1), DMDALocalIndex3D(info, i + 1, j + 0, k + 1), DMDALocalIndex3D(info, i + 1, j + 1, k + 1), DMDALocalIndex3D(info, i + 0, j + 1, k + 1)},
11659371c9d4SSatish Balay                          col2blocked[PRMNODE_SIZE * 4] = {DMDALocalIndex2D(info, i + 0, j + 0), DMDALocalIndex2D(info, i + 1, j + 0), DMDALocalIndex2D(info, i + 1, j + 1), DMDALocalIndex2D(info, i + 0, j + 1)};
1166c4762a1bSJed Brown #if !defined COMPUTE_LOWER_TRIANGULAR /* fill in lower-triangular part, this is really cheap compared to computing the entries */
1167c4762a1bSJed Brown           for (l = 0; l < 8; l++) {
1168c4762a1bSJed Brown             for (ll = l + 1; ll < 8; ll++) {
1169c4762a1bSJed Brown               Ke[ll * 2 + 0][l * 2 + 0] = Ke[l * 2 + 0][ll * 2 + 0];
1170c4762a1bSJed Brown               Ke[ll * 2 + 1][l * 2 + 0] = Ke[l * 2 + 0][ll * 2 + 1];
1171c4762a1bSJed Brown               Ke[ll * 2 + 0][l * 2 + 1] = Ke[l * 2 + 1][ll * 2 + 0];
1172c4762a1bSJed Brown               Ke[ll * 2 + 1][l * 2 + 1] = Ke[l * 2 + 1][ll * 2 + 1];
1173c4762a1bSJed Brown             }
1174c4762a1bSJed Brown           }
1175c4762a1bSJed Brown #endif
11769566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedLocal(B, 8, rc3blocked, 8, rc3blocked, &Ke[0][0], ADD_VALUES)); /* velocity-velocity coupling can use blocked insertion */
1177c4762a1bSJed Brown           {                                                                                            /* The off-diagonal part cannot (yet) */
1178c4762a1bSJed Brown             PetscInt row3scalar[NODE_SIZE * 8], col2scalar[PRMNODE_SIZE * 4];
11799371c9d4SSatish Balay             for (l = 0; l < 8; l++)
11809371c9d4SSatish Balay               for (ll = 0; ll < NODE_SIZE; ll++) row3scalar[l * NODE_SIZE + ll] = rc3blocked[l] * NODE_SIZE + ll;
11819371c9d4SSatish Balay             for (l = 0; l < 4; l++)
11829371c9d4SSatish Balay               for (ll = 0; ll < PRMNODE_SIZE; ll++) col2scalar[l * PRMNODE_SIZE + ll] = col2blocked[l] * PRMNODE_SIZE + ll;
11839566063dSJacob Faibussowitsch             PetscCall(MatSetValuesLocal(Bcpl, 8 * NODE_SIZE, row3scalar, 4 * PRMNODE_SIZE, col2scalar, &Kcpl[0][0], ADD_VALUES));
1184c4762a1bSJed Brown           }
1185c4762a1bSJed Brown         }
1186c4762a1bSJed Brown       }
1187c4762a1bSJed Brown     }
1188c4762a1bSJed Brown   }
1189c4762a1bSJed Brown   PetscFunctionReturn(0);
1190c4762a1bSJed Brown }
1191c4762a1bSJed Brown 
11929371c9d4SSatish Balay static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info, const Node ***x3, const PrmNode **x2, const PrmNode **xdot2, PetscReal a, Mat B22, Mat B21, THI thi) {
1193c4762a1bSJed Brown   PetscInt xs, ys, xm, ym, zm, i, j, k;
1194c4762a1bSJed Brown 
1195c4762a1bSJed Brown   PetscFunctionBeginUser;
1196c4762a1bSJed Brown   xs = info->zs;
1197c4762a1bSJed Brown   ys = info->ys;
1198c4762a1bSJed Brown   xm = info->zm;
1199c4762a1bSJed Brown   ym = info->ym;
1200c4762a1bSJed Brown   zm = info->xm;
1201c4762a1bSJed Brown 
12023c633725SBarry Smith   PetscCheck(zm <= 1024, ((PetscObject)info->da)->comm, PETSC_ERR_SUP, "Need to allocate more space");
1203c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
1204c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
1205c4762a1bSJed Brown       { /* Self-coupling */
1206c4762a1bSJed Brown         const PetscInt    row[]  = {DMDALocalIndex2D(info, i, j)};
1207c4762a1bSJed Brown         const PetscInt    col[]  = {DMDALocalIndex2D(info, i, j)};
12089371c9d4SSatish Balay         const PetscScalar vals[] = {a, 0, 0, 0, a, 0, 0, 0, a};
12099566063dSJacob Faibussowitsch         PetscCall(MatSetValuesBlockedLocal(B22, 1, row, 1, col, vals, INSERT_VALUES));
1210c4762a1bSJed Brown       }
1211c4762a1bSJed Brown       for (k = 0; k < zm; k++) { /* Coupling to velocity problem */
1212c4762a1bSJed Brown         /* Use a cheaper quadrature than for residual evaluation, because it is much sparser */
1213c4762a1bSJed Brown         const PetscInt    row[]  = {FieldIndex(PrmNode, DMDALocalIndex2D(info, i, j), h)};
12149371c9d4SSatish Balay         const PetscInt    cols[] = {FieldIndex(Node, DMDALocalIndex3D(info, i - 1, j, k), u), FieldIndex(Node, DMDALocalIndex3D(info, i, j, k), u), FieldIndex(Node, DMDALocalIndex3D(info, i + 1, j, k), u),
12159371c9d4SSatish Balay                                     FieldIndex(Node, DMDALocalIndex3D(info, i, j - 1, k), v), FieldIndex(Node, DMDALocalIndex3D(info, i, j, k), v), FieldIndex(Node, DMDALocalIndex3D(info, i, j + 1, k), v)};
12169371c9d4SSatish Balay         const PetscScalar w = (k && k < zm - 1) ? 0.5 : 0.25, hW = w * (x2[i - 1][j].h + x2[i][j].h) / (zm - 1.), hE = w * (x2[i][j].h + x2[i + 1][j].h) / (zm - 1.), hS = w * (x2[i][j - 1].h + x2[i][j].h) / (zm - 1.),
1217c4762a1bSJed Brown                           hN = w * (x2[i][j].h + x2[i][j + 1].h) / (zm - 1.);
12189371c9d4SSatish Balay         PetscScalar *vals, vals_upwind[] = {((PetscRealPart(x3[i][j][k].u) > 0) ? -hW : 0), ((PetscRealPart(x3[i][j][k].u) > 0) ? +hE : -hW), ((PetscRealPart(x3[i][j][k].u) > 0) ? 0 : +hE),
12199371c9d4SSatish Balay                                             ((PetscRealPart(x3[i][j][k].v) > 0) ? -hS : 0), ((PetscRealPart(x3[i][j][k].v) > 0) ? +hN : -hS), ((PetscRealPart(x3[i][j][k].v) > 0) ? 0 : +hN)},
12209371c9d4SSatish Balay                            vals_centered[] = {-0.5 * hW, 0.5 * (-hW + hE), 0.5 * hE, -0.5 * hS, 0.5 * (-hS + hN), 0.5 * hN};
1221c4762a1bSJed Brown         vals                               = 1 ? vals_upwind : vals_centered;
1222c4762a1bSJed Brown         if (k == 0) {
1223c4762a1bSJed Brown           Node derate;
1224c4762a1bSJed Brown           THIErosion(thi, &x3[i][j][0], NULL, &derate);
1225c4762a1bSJed Brown           vals[1] -= derate.u;
1226c4762a1bSJed Brown           vals[4] -= derate.v;
1227c4762a1bSJed Brown         }
12289566063dSJacob Faibussowitsch         PetscCall(MatSetValuesLocal(B21, 1, row, 6, cols, vals, INSERT_VALUES));
1229c4762a1bSJed Brown       }
1230c4762a1bSJed Brown     }
1231c4762a1bSJed Brown   }
1232c4762a1bSJed Brown   PetscFunctionReturn(0);
1233c4762a1bSJed Brown }
1234c4762a1bSJed Brown 
12359371c9d4SSatish Balay static PetscErrorCode THIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
1236c4762a1bSJed Brown   THI             thi = (THI)ctx;
1237c4762a1bSJed Brown   DM              pack, da3, da2;
1238c4762a1bSJed Brown   Vec             X3, X2, Xdot2;
1239c4762a1bSJed Brown   Mat             B11, B12, B21, B22;
1240c4762a1bSJed Brown   DMDALocalInfo   info3;
1241c4762a1bSJed Brown   IS             *isloc;
1242c4762a1bSJed Brown   const Node   ***x3;
1243c4762a1bSJed Brown   const PrmNode **x2, **xdot2;
1244c4762a1bSJed Brown 
1245c4762a1bSJed Brown   PetscFunctionBeginUser;
12469566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &pack));
12479566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntries(pack, &da3, &da2));
12489566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da3, &info3));
12499566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(pack, &X3, &X2));
12509566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(pack, NULL, &Xdot2));
12519566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(pack, X, X3, X2));
12529566063dSJacob Faibussowitsch   PetscCall(THIFixGhosts(thi, da3, da2, X3, X2));
12539566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(pack, Xdot, NULL, Xdot2));
1254c4762a1bSJed Brown 
12559566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
1256c4762a1bSJed Brown 
12579566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalISs(pack, &isloc));
12589566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(B, isloc[0], isloc[0], &B11));
12599566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(B, isloc[0], isloc[1], &B12));
12609566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(B, isloc[1], isloc[0], &B21));
12619566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(B, isloc[1], isloc[1], &B22));
1262c4762a1bSJed Brown 
12639566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da3, X3, &x3));
12649566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da2, X2, &x2));
12659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da2, Xdot2, &xdot2));
1266c4762a1bSJed Brown 
12679566063dSJacob Faibussowitsch   PetscCall(THIJacobianLocal_Momentum(&info3, x3, x2, B11, B12, thi));
1268c4762a1bSJed Brown 
1269c4762a1bSJed Brown   /* Need to switch from ADD_VALUES to INSERT_VALUES */
12709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FLUSH_ASSEMBLY));
12719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FLUSH_ASSEMBLY));
1272c4762a1bSJed Brown 
12739566063dSJacob Faibussowitsch   PetscCall(THIJacobianLocal_2D(&info3, x3, x2, xdot2, a, B22, B21, thi));
1274c4762a1bSJed Brown 
12759566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da3, X3, &x3));
12769566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da2, X2, &x2));
12779566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da2, Xdot2, &xdot2));
1278c4762a1bSJed Brown 
12799566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(B, isloc[0], isloc[0], &B11));
12809566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(B, isloc[0], isloc[1], &B12));
12819566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(B, isloc[1], isloc[0], &B21));
12829566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(B, isloc[1], isloc[1], &B22));
12839566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isloc[0]));
12849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isloc[1]));
12859566063dSJacob Faibussowitsch   PetscCall(PetscFree(isloc));
1286c4762a1bSJed Brown 
12879566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(pack, &X3, &X2));
12889566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(pack, 0, &Xdot2));
1289c4762a1bSJed Brown 
12909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
12919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1292c4762a1bSJed Brown   if (A != B) {
12939566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
12949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1295c4762a1bSJed Brown   }
12969566063dSJacob Faibussowitsch   if (thi->verbose) PetscCall(THIMatrixStatistics(thi, B, PETSC_VIEWER_STDOUT_WORLD));
1297c4762a1bSJed Brown   PetscFunctionReturn(0);
1298c4762a1bSJed Brown }
1299c4762a1bSJed Brown 
1300c4762a1bSJed Brown /* VTK's XML formats are so brain-dead that they can't handle multiple grids in the same file.  Since the communication
1301c4762a1bSJed Brown  * can be shared between the two grids, we write two files at once, one for velocity (living on a 3D grid defined by
1302c4762a1bSJed Brown  * h=thickness and b=bed) and another for all properties living on the 2D grid.
1303c4762a1bSJed Brown  */
13049371c9d4SSatish Balay static PetscErrorCode THIDAVecView_VTK_XML(THI thi, DM pack, Vec X, const char filename[], const char filename2[]) {
1305c4762a1bSJed Brown   const PetscInt dof = NODE_SIZE, dof2 = PRMNODE_SIZE;
1306c4762a1bSJed Brown   Units          units = thi->units;
1307c4762a1bSJed Brown   MPI_Comm       comm;
1308c4762a1bSJed Brown   PetscViewer    viewer3, viewer2;
1309c4762a1bSJed Brown   PetscMPIInt    rank, size, tag, nn, nmax, nn2, nmax2;
1310c4762a1bSJed Brown   PetscInt       mx, my, mz, r, range[6];
1311c4762a1bSJed Brown   PetscScalar   *x, *x2;
1312c4762a1bSJed Brown   DM             da3, da2;
1313c4762a1bSJed Brown   Vec            X3, X2;
1314c4762a1bSJed Brown 
1315c4762a1bSJed Brown   PetscFunctionBeginUser;
13169566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)thi, &comm));
13179566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntries(pack, &da3, &da2));
13189566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccess(pack, X, &X3, &X2));
13199566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da3, 0, &mz, &my, &mx, 0, 0, 0, 0, 0, 0, 0, 0, 0));
13209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
13219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
13229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIOpen(comm, filename, &viewer3));
13239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIOpen(comm, filename2, &viewer2));
13249566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer3, "<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n"));
13259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer2, "<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n"));
132663a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer3, "  <StructuredGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n", 0, mz - 1, 0, my - 1, 0, mx - 1));
132763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer2, "  <StructuredGrid WholeExtent=\"%d %d %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n", 0, 0, 0, my - 1, 0, mx - 1));
1328c4762a1bSJed Brown 
13299566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da3, range, range + 1, range + 2, range + 3, range + 4, range + 5));
13309566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(range[3] * range[4] * range[5] * dof, &nn));
13319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&nn, &nmax, 1, MPI_INT, MPI_MAX, 0, comm));
13329566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(range[4] * range[5] * dof2, &nn2));
13339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&nn2, &nmax2, 1, MPI_INT, MPI_MAX, 0, comm));
1334c4762a1bSJed Brown   tag = ((PetscObject)viewer3)->tag;
13359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X3, (const PetscScalar **)&x));
13369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X2, (const PetscScalar **)&x2));
1337dd400576SPatrick Sanan   if (rank == 0) {
1338c4762a1bSJed Brown     PetscScalar *array, *array2;
13399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nmax, &array, nmax2, &array2));
1340c4762a1bSJed Brown     for (r = 0; r < size; r++) {
1341c4762a1bSJed Brown       PetscInt i, j, k, f, xs, xm, ys, ym, zs, zm;
1342c4762a1bSJed Brown       Node    *y3;
1343c4762a1bSJed Brown       PetscScalar(*y2)[PRMNODE_SIZE];
1344c4762a1bSJed Brown       MPI_Status status;
1345c4762a1bSJed Brown 
1346*48a46eb9SPierre Jolivet       if (r) PetscCallMPI(MPI_Recv(range, 6, MPIU_INT, r, tag, comm, MPI_STATUS_IGNORE));
13479371c9d4SSatish Balay       zs = range[0];
13489371c9d4SSatish Balay       ys = range[1];
13499371c9d4SSatish Balay       xs = range[2];
13509371c9d4SSatish Balay       zm = range[3];
13519371c9d4SSatish Balay       ym = range[4];
13529371c9d4SSatish Balay       xm = range[5];
13533c633725SBarry Smith       PetscCheck(xm * ym * zm * dof <= nmax, PETSC_COMM_SELF, PETSC_ERR_PLIB, "should not happen");
1354c4762a1bSJed Brown       if (r) {
13559566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(array, nmax, MPIU_SCALAR, r, tag, comm, &status));
13569566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn));
13573c633725SBarry Smith         PetscCheck(nn == xm * ym * zm * dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "corrupt da3 send");
1358c4762a1bSJed Brown         y3 = (Node *)array;
13599566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(array2, nmax2, MPIU_SCALAR, r, tag, comm, &status));
13609566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn2));
13613c633725SBarry Smith         PetscCheck(nn2 == xm * ym * dof2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "corrupt da2 send");
1362c4762a1bSJed Brown         y2 = (PetscScalar(*)[PRMNODE_SIZE])array2;
1363c4762a1bSJed Brown       } else {
1364c4762a1bSJed Brown         y3 = (Node *)x;
1365c4762a1bSJed Brown         y2 = (PetscScalar(*)[PRMNODE_SIZE])x2;
1366c4762a1bSJed Brown       }
136763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer3, "    <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n", zs, zs + zm - 1, ys, ys + ym - 1, xs, xs + xm - 1));
136863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer2, "    <Piece Extent=\"%d %d %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n", 0, 0, ys, ys + ym - 1, xs, xs + xm - 1));
1369c4762a1bSJed Brown 
13709566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer3, "      <Points>\n"));
13719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer2, "      <Points>\n"));
13729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer3, "        <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n"));
13739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer2, "        <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n"));
1374c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
1375c4762a1bSJed Brown         for (j = ys; j < ys + ym; j++) {
13769371c9d4SSatish Balay           PetscReal xx = thi->Lx * i / mx, yy = thi->Ly * j / my, b = PetscRealPart(y2[i * ym + j][FieldOffset(PrmNode, b)]), h = PetscRealPart(y2[i * ym + j][FieldOffset(PrmNode, h)]);
1377c4762a1bSJed Brown           for (k = zs; k < zs + zm; k++) {
1378c4762a1bSJed Brown             PetscReal zz = b + h * k / (mz - 1.);
137963a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer3, "%f %f %f\n", (double)xx, (double)yy, (double)zz));
1380c4762a1bSJed Brown           }
138163a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer2, "%f %f %f\n", (double)xx, (double)yy, (double)0.0));
1382c4762a1bSJed Brown         }
1383c4762a1bSJed Brown       }
13849566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer3, "        </DataArray>\n"));
13859566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer2, "        </DataArray>\n"));
13869566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer3, "      </Points>\n"));
13879566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer2, "      </Points>\n"));
1388c4762a1bSJed Brown 
1389c4762a1bSJed Brown       { /* Velocity and rank (3D) */
13909566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer3, "      <PointData>\n"));
13919566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer3, "        <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n"));
1392*48a46eb9SPierre Jolivet         for (i = 0; i < nn / dof; i++) PetscCall(PetscViewerASCIIPrintf(viewer3, "%f %f %f\n", (double)(PetscRealPart(y3[i].u) * units->year / units->meter), (double)(PetscRealPart(y3[i].v) * units->year / units->meter), 0.0));
13939566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer3, "        </DataArray>\n"));
1394c4762a1bSJed Brown 
13959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer3, "        <DataArray type=\"Int32\" Name=\"rank\" NumberOfComponents=\"1\" format=\"ascii\">\n"));
1396*48a46eb9SPierre Jolivet         for (i = 0; i < nn; i += dof) PetscCall(PetscViewerASCIIPrintf(viewer3, "%" PetscInt_FMT "\n", r));
13979566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer3, "        </DataArray>\n"));
13989566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer3, "      </PointData>\n"));
1399c4762a1bSJed Brown       }
1400c4762a1bSJed Brown 
1401c4762a1bSJed Brown       { /* 2D */
14029566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer2, "      <PointData>\n"));
1403c4762a1bSJed Brown         for (f = 0; f < PRMNODE_SIZE; f++) {
1404c4762a1bSJed Brown           const char *fieldname;
14059566063dSJacob Faibussowitsch           PetscCall(DMDAGetFieldName(da2, f, &fieldname));
14069566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer2, "        <DataArray type=\"Float32\" Name=\"%s\" format=\"ascii\">\n", fieldname));
1407*48a46eb9SPierre Jolivet           for (i = 0; i < nn2 / PRMNODE_SIZE; i++) PetscCall(PetscViewerASCIIPrintf(viewer2, "%g\n", (double)y2[i][f]));
14089566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer2, "        </DataArray>\n"));
1409c4762a1bSJed Brown         }
14109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer2, "      </PointData>\n"));
1411c4762a1bSJed Brown       }
1412c4762a1bSJed Brown 
14139566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer3, "    </Piece>\n"));
14149566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer2, "    </Piece>\n"));
1415c4762a1bSJed Brown     }
14169566063dSJacob Faibussowitsch     PetscCall(PetscFree2(array, array2));
1417c4762a1bSJed Brown   } else {
14189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(range, 6, MPIU_INT, 0, tag, comm));
14199566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(x, nn, MPIU_SCALAR, 0, tag, comm));
14209566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(x2, nn2, MPIU_SCALAR, 0, tag, comm));
1421c4762a1bSJed Brown   }
14229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X3, (const PetscScalar **)&x));
14239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X2, (const PetscScalar **)&x2));
14249566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer3, "  </StructuredGrid>\n"));
14259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer2, "  </StructuredGrid>\n"));
1426c4762a1bSJed Brown 
14279566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccess(pack, X, &X3, &X2));
14289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer3, "</VTKFile>\n"));
14299566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer2, "</VTKFile>\n"));
14309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer3));
14319566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer2));
1432c4762a1bSJed Brown   PetscFunctionReturn(0);
1433c4762a1bSJed Brown }
1434c4762a1bSJed Brown 
14359371c9d4SSatish Balay static PetscErrorCode THITSMonitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
1436c4762a1bSJed Brown   THI  thi = (THI)ctx;
1437c4762a1bSJed Brown   DM   pack;
1438c4762a1bSJed Brown   char filename3[PETSC_MAX_PATH_LEN], filename2[PETSC_MAX_PATH_LEN];
1439c4762a1bSJed Brown 
1440c4762a1bSJed Brown   PetscFunctionBeginUser;
1441c4762a1bSJed Brown   if (step < 0) PetscFunctionReturn(0); /* negative one is used to indicate an interpolated solution */
144263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%3" PetscInt_FMT ": t=%g\n", step, (double)t));
1443c4762a1bSJed Brown   if (thi->monitor_interval && step % thi->monitor_interval) PetscFunctionReturn(0);
14449566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &pack));
144563a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename3, sizeof(filename3), "%s-3d-%03" PetscInt_FMT ".vts", thi->monitor_basename, step));
144663a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename2, sizeof(filename2), "%s-2d-%03" PetscInt_FMT ".vts", thi->monitor_basename, step));
14479566063dSJacob Faibussowitsch   PetscCall(THIDAVecView_VTK_XML(thi, pack, X, filename3, filename2));
1448c4762a1bSJed Brown   PetscFunctionReturn(0);
1449c4762a1bSJed Brown }
1450c4762a1bSJed Brown 
14519371c9d4SSatish Balay static PetscErrorCode THICreateDM3d(THI thi, DM *dm3d) {
1452c4762a1bSJed Brown   MPI_Comm comm;
1453c4762a1bSJed Brown   PetscInt M = 3, N = 3, P = 2;
1454c4762a1bSJed Brown   DM       da;
1455c4762a1bSJed Brown 
1456c4762a1bSJed Brown   PetscFunctionBeginUser;
14579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)thi, &comm));
1458d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Grid resolution options", "");
1459c4762a1bSJed Brown   {
14609566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-M", "Number of elements in x-direction on coarse level", "", M, &M, NULL));
1461c4762a1bSJed Brown     N = M;
14629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-N", "Number of elements in y-direction on coarse level (if different from M)", "", N, &N, NULL));
14639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-P", "Number of elements in z-direction on coarse level", "", P, &P, NULL));
1464c4762a1bSJed Brown   }
1465d0609cedSBarry Smith   PetscOptionsEnd();
14669566063dSJacob Faibussowitsch   PetscCall(DMDACreate3d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, DMDA_STENCIL_BOX, P, N, M, 1, PETSC_DETERMINE, PETSC_DETERMINE, sizeof(Node) / sizeof(PetscScalar), 1, 0, 0, 0, &da));
14679566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
14689566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
14699566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "x-velocity"));
14709566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "y-velocity"));
1471c4762a1bSJed Brown   *dm3d = da;
1472c4762a1bSJed Brown   PetscFunctionReturn(0);
1473c4762a1bSJed Brown }
1474c4762a1bSJed Brown 
14759371c9d4SSatish Balay int main(int argc, char *argv[]) {
1476c4762a1bSJed Brown   MPI_Comm  comm;
1477c4762a1bSJed Brown   DM        pack, da3, da2;
1478c4762a1bSJed Brown   TS        ts;
1479c4762a1bSJed Brown   THI       thi;
1480c4762a1bSJed Brown   Vec       X;
1481c4762a1bSJed Brown   Mat       B;
1482c4762a1bSJed Brown   PetscInt  i, steps;
1483c4762a1bSJed Brown   PetscReal ftime;
1484c4762a1bSJed Brown 
1485327415f7SBarry Smith   PetscFunctionBeginUser;
14869566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
1487c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1488c4762a1bSJed Brown 
14899566063dSJacob Faibussowitsch   PetscCall(THICreate(comm, &thi));
14909566063dSJacob Faibussowitsch   PetscCall(THICreateDM3d(thi, &da3));
1491c4762a1bSJed Brown   {
1492c4762a1bSJed Brown     PetscInt        Mx, My, mx, my, s;
1493c4762a1bSJed Brown     DMDAStencilType st;
14949566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da3, 0, 0, &My, &Mx, 0, &my, &mx, 0, &s, 0, 0, 0, &st));
14959566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)thi), DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, st, My, Mx, my, mx, sizeof(PrmNode) / sizeof(PetscScalar), s, 0, 0, &da2));
14969566063dSJacob Faibussowitsch     PetscCall(DMSetUp(da2));
1497c4762a1bSJed Brown   }
1498c4762a1bSJed Brown 
14999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)da3, "3D_Velocity"));
15009566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(da3, "f3d_"));
15019566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da3, 0, "u"));
15029566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da3, 1, "v"));
15039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)da2, "2D_Fields"));
15049566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(da2, "f2d_"));
15059566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da2, 0, "b"));
15069566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da2, 1, "h"));
15079566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da2, 2, "beta2"));
15089566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, &pack));
15099566063dSJacob Faibussowitsch   PetscCall(DMCompositeAddDM(pack, da3));
15109566063dSJacob Faibussowitsch   PetscCall(DMCompositeAddDM(pack, da2));
15119566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da3));
15129566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da2));
15139566063dSJacob Faibussowitsch   PetscCall(DMSetUp(pack));
15149566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pack, &B));
15159566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
15169566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(B, "thi_"));
1517c4762a1bSJed Brown 
1518c4762a1bSJed Brown   for (i = 0; i < thi->nlevels; i++) {
1519c4762a1bSJed Brown     PetscReal Lx = thi->Lx / thi->units->meter, Ly = thi->Ly / thi->units->meter, Lz = thi->Lz / thi->units->meter;
1520c4762a1bSJed Brown     PetscInt  Mx, My, Mz;
15219566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntries(pack, &da3, &da2));
15229566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da3, 0, &Mz, &My, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0));
152363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi), "Level %" PetscInt_FMT " domain size (m) %8.2g x %8.2g x %8.2g, num elements %3d x %3d x %3d (%8d), size (m) %g x %g x %g\n", i, Lx, Ly, Lz, Mx, My, Mz, Mx * My * Mz, Lx / Mx, Ly / My, 1000. / (Mz - 1)));
1524c4762a1bSJed Brown   }
1525c4762a1bSJed Brown 
15269566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(pack, &X));
15279566063dSJacob Faibussowitsch   PetscCall(THIInitial(thi, pack, X));
1528c4762a1bSJed Brown 
15299566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
15309566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, pack));
15319566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
15329566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, THITSMonitor, thi, NULL));
15339566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSTHETA));
15349566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, THIFunction, thi));
15359566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, B, B, THIJacobian, thi));
15369566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10.0));
15379566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
15389566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
15399566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1e-3));
15409566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1541c4762a1bSJed Brown 
15429566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
15439566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
15449566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
154563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Steps %" PetscInt_FMT "  final time %g\n", steps, (double)ftime));
1546c4762a1bSJed Brown 
15479566063dSJacob Faibussowitsch   if (0) PetscCall(THISolveStatistics(thi, ts, 0, "Full"));
1548c4762a1bSJed Brown 
1549c4762a1bSJed Brown   {
1550c4762a1bSJed Brown     PetscBool flg;
1551c4762a1bSJed Brown     char      filename[PETSC_MAX_PATH_LEN] = "";
15529566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(NULL, NULL, "-o", filename, sizeof(filename), &flg));
15531baa6e33SBarry Smith     if (flg) PetscCall(THIDAVecView_VTK_XML(thi, pack, X, filename, NULL));
1554c4762a1bSJed Brown   }
1555c4762a1bSJed Brown 
15569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
15579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
15589566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pack));
15599566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
15609566063dSJacob Faibussowitsch   PetscCall(THIDestroy(&thi));
15619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1562b122ec5aSJacob Faibussowitsch   return 0;
1563c4762a1bSJed Brown }
1564