xref: /petsc/src/ts/tutorials/ex14.c (revision 2f613bf53f46f9356e00a2ca2bd69453be72fc31)
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 */
58c4762a1bSJed Brown #define USE_SSE2_KERNELS (!defined NO_SSE2                              \
59c4762a1bSJed Brown                           && !defined PETSC_USE_COMPLEX                 \
60c4762a1bSJed Brown                           && !defined PETSC_USE_REAL_SINGLE           \
61c4762a1bSJed Brown                           && defined __SSE2__)
62c4762a1bSJed Brown 
63c4762a1bSJed Brown #if !defined __STDC_VERSION__ || __STDC_VERSION__ < 199901L
64c4762a1bSJed Brown #  if defined __cplusplus       /* C++ restrict is nonstandard and compilers have inconsistent rules about where it can be used */
65c4762a1bSJed Brown #    define restrict
66c4762a1bSJed Brown #  else
67c4762a1bSJed Brown #    define restrict PETSC_RESTRICT
68c4762a1bSJed Brown #  endif
69c4762a1bSJed Brown #endif
70c4762a1bSJed Brown 
71c4762a1bSJed Brown static PetscClassId THI_CLASSID;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown typedef enum {QUAD_GAUSS,QUAD_LOBATTO} 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 */
83c4762a1bSJed Brown static const PetscReal HexQInterp_Lobatto[8][8] = {{H,0,0,0,L,0,0,0},
84c4762a1bSJed Brown                                                    {0,H,0,0,0,L,0,0},
85c4762a1bSJed Brown                                                    {0,0,H,0,0,0,L,0},
86c4762a1bSJed Brown                                                    {0,0,0,H,0,0,0,L},
87c4762a1bSJed Brown                                                    {L,0,0,0,H,0,0,0},
88c4762a1bSJed Brown                                                    {0,L,0,0,0,H,0,0},
89c4762a1bSJed Brown                                                    {0,0,L,0,0,0,H,0},
90c4762a1bSJed Brown                                                    {0,0,0,L,0,0,0,H}};
91c4762a1bSJed Brown static const PetscReal HexQDeriv_Lobatto[8][8][3] = {
92c4762a1bSJed 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}  },
93c4762a1bSJed 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}    },
94c4762a1bSJed 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}  },
95c4762a1bSJed 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}},
96c4762a1bSJed 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}  },
97c4762a1bSJed 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}    },
98c4762a1bSJed 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}  },
99c4762a1bSJed Brown   {{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}}};
100c4762a1bSJed Brown /* Stanndard Gauss */
101c4762a1bSJed Brown static const PetscReal HexQInterp_Gauss[8][8] = {{H*H*H,L*H*H,L*L*H,H*L*H, H*H*L,L*H*L,L*L*L,H*L*L},
102c4762a1bSJed 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},
103c4762a1bSJed 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},
104c4762a1bSJed 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},
105c4762a1bSJed 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},
106c4762a1bSJed 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},
107c4762a1bSJed 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},
108c4762a1bSJed Brown                                                  {H*L*L,L*L*L,L*H*L,H*H*L, H*L*H,L*L*H,L*H*H,H*H*H}};
109c4762a1bSJed Brown static const PetscReal HexQDeriv_Gauss[8][8][3] = {
110c4762a1bSJed 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}},
111c4762a1bSJed 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}},
112c4762a1bSJed 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}},
113c4762a1bSJed 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}},
114c4762a1bSJed 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}},
115c4762a1bSJed 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}},
116c4762a1bSJed 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}},
117c4762a1bSJed Brown   {{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}}};
118c4762a1bSJed Brown static const PetscReal (*HexQInterp)[8],(*HexQDeriv)[8][3];
119c4762a1bSJed Brown /* Standard 2x2 Gauss quadrature for the bottom layer. */
120c4762a1bSJed Brown static const PetscReal QuadQInterp[4][4] = {{H*H,L*H,L*L,H*L},
121c4762a1bSJed Brown                                             {L*H,H*H,H*L,L*L},
122c4762a1bSJed Brown                                             {L*L,H*L,H*H,L*H},
123c4762a1bSJed Brown                                             {H*L,L*L,L*H,H*H}};
124c4762a1bSJed Brown static const PetscReal QuadQDeriv[4][4][2] = {
125c4762a1bSJed Brown   {{M*H,M*H},{P*H,M*L},{P*L,P*L},{M*L,P*H}},
126c4762a1bSJed Brown   {{M*H,M*L},{P*H,M*H},{P*L,P*H},{M*L,P*L}},
127c4762a1bSJed Brown   {{M*L,M*L},{P*L,M*H},{P*H,P*H},{M*H,P*L}},
128c4762a1bSJed Brown   {{M*L,M*H},{P*L,M*L},{P*H,P*L},{M*H,P*H}}};
129c4762a1bSJed Brown #undef G
130c4762a1bSJed Brown #undef H
131c4762a1bSJed Brown #undef L
132c4762a1bSJed Brown #undef M
133c4762a1bSJed Brown #undef P
134c4762a1bSJed Brown 
135c4762a1bSJed Brown #define HexExtract(x,i,j,k,n) do {              \
136c4762a1bSJed Brown     (n)[0] = (x)[i][j][k];                      \
137c4762a1bSJed Brown     (n)[1] = (x)[i+1][j][k];                    \
138c4762a1bSJed Brown     (n)[2] = (x)[i+1][j+1][k];                  \
139c4762a1bSJed Brown     (n)[3] = (x)[i][j+1][k];                    \
140c4762a1bSJed Brown     (n)[4] = (x)[i][j][k+1];                    \
141c4762a1bSJed Brown     (n)[5] = (x)[i+1][j][k+1];                  \
142c4762a1bSJed Brown     (n)[6] = (x)[i+1][j+1][k+1];                \
143c4762a1bSJed Brown     (n)[7] = (x)[i][j+1][k+1];                  \
144c4762a1bSJed Brown   } while (0)
145c4762a1bSJed Brown 
146c4762a1bSJed Brown #define HexExtractRef(x,i,j,k,n) do {           \
147c4762a1bSJed Brown     (n)[0] = &(x)[i][j][k];                     \
148c4762a1bSJed Brown     (n)[1] = &(x)[i+1][j][k];                   \
149c4762a1bSJed Brown     (n)[2] = &(x)[i+1][j+1][k];                 \
150c4762a1bSJed Brown     (n)[3] = &(x)[i][j+1][k];                   \
151c4762a1bSJed Brown     (n)[4] = &(x)[i][j][k+1];                   \
152c4762a1bSJed Brown     (n)[5] = &(x)[i+1][j][k+1];                 \
153c4762a1bSJed Brown     (n)[6] = &(x)[i+1][j+1][k+1];               \
154c4762a1bSJed Brown     (n)[7] = &(x)[i][j+1][k+1];                 \
155c4762a1bSJed Brown   } while (0)
156c4762a1bSJed Brown 
157c4762a1bSJed Brown #define QuadExtract(x,i,j,n) do {               \
158c4762a1bSJed Brown     (n)[0] = (x)[i][j];                         \
159c4762a1bSJed Brown     (n)[1] = (x)[i+1][j];                       \
160c4762a1bSJed Brown     (n)[2] = (x)[i+1][j+1];                     \
161c4762a1bSJed Brown     (n)[3] = (x)[i][j+1];                       \
162c4762a1bSJed Brown   } while (0)
163c4762a1bSJed Brown 
164c4762a1bSJed Brown static PetscScalar Sqr(PetscScalar a) {return a*a;}
165c4762a1bSJed Brown 
166c4762a1bSJed Brown static void HexGrad(const PetscReal dphi[][3],const PetscReal zn[],PetscReal dz[])
167c4762a1bSJed Brown {
168c4762a1bSJed Brown   PetscInt i;
169c4762a1bSJed Brown   dz[0] = dz[1] = dz[2] = 0;
170c4762a1bSJed Brown   for (i=0; i<8; i++) {
171c4762a1bSJed Brown     dz[0] += dphi[i][0] * zn[i];
172c4762a1bSJed Brown     dz[1] += dphi[i][1] * zn[i];
173c4762a1bSJed Brown     dz[2] += dphi[i][2] * zn[i];
174c4762a1bSJed Brown   }
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown static void HexComputeGeometry(PetscInt q,PetscReal hx,PetscReal hy,const PetscReal dz[restrict],PetscReal phi[restrict],PetscReal dphi[restrict][3],PetscReal *restrict jw)
178c4762a1bSJed Brown {
179c4762a1bSJed Brown   const PetscReal
180c4762a1bSJed Brown     jac[3][3] = {{hx/2,0,0}, {0,hy/2,0}, {dz[0],dz[1],dz[2]}}
181c4762a1bSJed Brown   ,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]}}
182c4762a1bSJed Brown   ,jdet = jac[0][0]*jac[1][1]*jac[2][2];
183c4762a1bSJed Brown   PetscInt i;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   for (i=0; i<8; i++) {
186c4762a1bSJed Brown     const PetscReal *dphir = HexQDeriv[q][i];
187c4762a1bSJed Brown     phi[i] = HexQInterp[q][i];
188c4762a1bSJed Brown     dphi[i][0] = dphir[0]*ijac[0][0] + dphir[1]*ijac[1][0] + dphir[2]*ijac[2][0];
189c4762a1bSJed Brown     dphi[i][1] = dphir[0]*ijac[0][1] + dphir[1]*ijac[1][1] + dphir[2]*ijac[2][1];
190c4762a1bSJed Brown     dphi[i][2] = dphir[0]*ijac[0][2] + dphir[1]*ijac[1][2] + dphir[2]*ijac[2][2];
191c4762a1bSJed Brown   }
192c4762a1bSJed Brown   *jw = 1.0 * jdet;
193c4762a1bSJed Brown }
194c4762a1bSJed Brown 
195c4762a1bSJed Brown typedef struct _p_THI   *THI;
196c4762a1bSJed Brown typedef struct _n_Units *Units;
197c4762a1bSJed Brown 
198c4762a1bSJed Brown typedef struct {
199c4762a1bSJed Brown   PetscScalar u,v;
200c4762a1bSJed Brown } Node;
201c4762a1bSJed Brown 
202c4762a1bSJed Brown typedef struct {
203c4762a1bSJed Brown   PetscScalar b;                /* bed */
204c4762a1bSJed Brown   PetscScalar h;                /* thickness */
205c4762a1bSJed Brown   PetscScalar beta2;            /* friction */
206c4762a1bSJed Brown } PrmNode;
207c4762a1bSJed Brown 
208c4762a1bSJed Brown #define FieldSize(ntype) ((PetscInt)(sizeof(ntype)/sizeof(PetscScalar)))
209c4762a1bSJed Brown #define FieldOffset(ntype,member) ((PetscInt)(offsetof(ntype,member)/sizeof(PetscScalar)))
210c4762a1bSJed Brown #define FieldIndex(ntype,i,member) ((PetscInt)((i)*FieldSize(ntype) + FieldOffset(ntype,member)))
211c4762a1bSJed Brown #define NODE_SIZE FieldSize(Node)
212c4762a1bSJed Brown #define PRMNODE_SIZE FieldSize(PrmNode)
213c4762a1bSJed Brown 
214c4762a1bSJed Brown typedef struct {
215c4762a1bSJed Brown   PetscReal min,max,cmin,cmax;
216c4762a1bSJed Brown } PRange;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown struct _p_THI {
219c4762a1bSJed Brown   PETSCHEADER(int);
220c4762a1bSJed Brown   void      (*initialize)(THI,PetscReal x,PetscReal y,PrmNode *p);
221c4762a1bSJed Brown   PetscInt  nlevels;
222c4762a1bSJed Brown   PetscInt  zlevels;
223c4762a1bSJed Brown   PetscReal Lx,Ly,Lz;           /* Model domain */
224c4762a1bSJed Brown   PetscReal alpha;              /* Bed angle */
225c4762a1bSJed Brown   Units     units;
226c4762a1bSJed Brown   PetscReal dirichlet_scale;
227c4762a1bSJed Brown   PetscReal ssa_friction_scale;
228c4762a1bSJed Brown   PetscReal inertia;
229c4762a1bSJed Brown   PRange    eta;
230c4762a1bSJed Brown   PRange    beta2;
231c4762a1bSJed Brown   struct {
232c4762a1bSJed Brown     PetscReal Bd2,eps,exponent,glen_n;
233c4762a1bSJed Brown   } viscosity;
234c4762a1bSJed Brown   struct {
235c4762a1bSJed Brown     PetscReal irefgam,eps2,exponent;
236c4762a1bSJed Brown   } friction;
237c4762a1bSJed Brown   struct {
238c4762a1bSJed Brown     PetscReal rate,exponent,refvel;
239c4762a1bSJed Brown   } erosion;
240c4762a1bSJed Brown   PetscReal rhog;
241c4762a1bSJed Brown   PetscBool no_slip;
242c4762a1bSJed Brown   PetscBool verbose;
243c4762a1bSJed Brown   char      *mattype;
244c4762a1bSJed Brown   char      *monitor_basename;
245c4762a1bSJed Brown   PetscInt  monitor_interval;
246c4762a1bSJed Brown };
247c4762a1bSJed Brown 
248c4762a1bSJed Brown struct _n_Units {
249c4762a1bSJed Brown   /* fundamental */
250c4762a1bSJed Brown   PetscReal meter;
251c4762a1bSJed Brown   PetscReal kilogram;
252c4762a1bSJed Brown   PetscReal second;
253c4762a1bSJed Brown   /* derived */
254c4762a1bSJed Brown   PetscReal Pascal;
255c4762a1bSJed Brown   PetscReal year;
256c4762a1bSJed Brown };
257c4762a1bSJed Brown 
258c4762a1bSJed Brown static void PrmHexGetZ(const PrmNode pn[],PetscInt k,PetscInt zm,PetscReal zn[])
259c4762a1bSJed Brown {
260c4762a1bSJed Brown   const PetscScalar zm1 = zm-1,
261c4762a1bSJed Brown     znl[8] = {pn[0].b + pn[0].h*(PetscScalar)k/zm1,
262c4762a1bSJed Brown               pn[1].b + pn[1].h*(PetscScalar)k/zm1,
263c4762a1bSJed Brown               pn[2].b + pn[2].h*(PetscScalar)k/zm1,
264c4762a1bSJed Brown               pn[3].b + pn[3].h*(PetscScalar)k/zm1,
265c4762a1bSJed Brown               pn[0].b + pn[0].h*(PetscScalar)(k+1)/zm1,
266c4762a1bSJed Brown               pn[1].b + pn[1].h*(PetscScalar)(k+1)/zm1,
267c4762a1bSJed Brown               pn[2].b + pn[2].h*(PetscScalar)(k+1)/zm1,
268c4762a1bSJed Brown               pn[3].b + pn[3].h*(PetscScalar)(k+1)/zm1};
269c4762a1bSJed Brown   PetscInt i;
270c4762a1bSJed Brown   for (i=0; i<8; i++) zn[i] = PetscRealPart(znl[i]);
271c4762a1bSJed Brown }
272c4762a1bSJed Brown 
273c4762a1bSJed Brown /* Compute a gradient of all the 2D fields at four quadrature points.  Output for [quadrature_point][direction].field_name */
274c4762a1bSJed Brown static PetscErrorCode QuadComputeGrad4(const PetscReal dphi[][4][2],PetscReal hx,PetscReal hy,const PrmNode pn[4],PrmNode dp[4][2])
275c4762a1bSJed Brown {
276c4762a1bSJed Brown   PetscErrorCode ierr;
277c4762a1bSJed Brown   PetscInt       q,i,f;
278c4762a1bSJed Brown   const PetscScalar (*restrict pg)[PRMNODE_SIZE] = (const PetscScalar(*)[PRMNODE_SIZE])pn; /* Get generic array pointers to the node */
279c4762a1bSJed Brown   PetscScalar (*restrict dpg)[2][PRMNODE_SIZE]   = (PetscScalar(*)[2][PRMNODE_SIZE])dp;
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   PetscFunctionBeginUser;
282c4762a1bSJed Brown   ierr = PetscArrayzero(dpg,4);CHKERRQ(ierr);
283c4762a1bSJed Brown   for (q=0; q<4; q++) {
284c4762a1bSJed Brown     for (i=0; i<4; i++) {
285c4762a1bSJed Brown       for (f=0; f<PRMNODE_SIZE; f++) {
286c4762a1bSJed Brown         dpg[q][0][f] += dphi[q][i][0]/hx * pg[i][f];
287c4762a1bSJed Brown         dpg[q][1][f] += dphi[q][i][1]/hy * pg[i][f];
288c4762a1bSJed Brown       }
289c4762a1bSJed Brown     }
290c4762a1bSJed Brown   }
291c4762a1bSJed Brown   PetscFunctionReturn(0);
292c4762a1bSJed Brown }
293c4762a1bSJed Brown 
294c4762a1bSJed Brown static inline PetscReal StaggeredMidpoint2D(PetscScalar a,PetscScalar b,PetscScalar c,PetscScalar d)
295c4762a1bSJed Brown {return 0.5*PetscRealPart(0.75*a + 0.75*b + 0.25*c + 0.25*d);}
296c4762a1bSJed Brown static inline PetscReal UpwindFlux1D(PetscReal u,PetscReal hL,PetscReal hR)
297c4762a1bSJed Brown {return (u > 0) ? hL*u : hR*u;}
298c4762a1bSJed Brown 
299c4762a1bSJed Brown #define UpwindFluxXW(x3,x2,h,i,j,k,dj) 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), \
300c4762a1bSJed Brown                                                     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))
301c4762a1bSJed Brown #define UpwindFluxXE(x3,x2,h,i,j,k,dj) 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), \
302c4762a1bSJed Brown                                                     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))
303c4762a1bSJed Brown #define UpwindFluxYS(x3,x2,h,i,j,k,di) 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), \
304c4762a1bSJed Brown                                                     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))
305c4762a1bSJed Brown #define UpwindFluxYN(x3,x2,h,i,j,k,di) 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), \
306c4762a1bSJed Brown                                                     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))
307c4762a1bSJed Brown 
308c4762a1bSJed Brown static void PrmNodeGetFaceMeasure(const PrmNode **p,PetscInt i,PetscInt j,PetscScalar h[])
309c4762a1bSJed Brown {
310c4762a1bSJed Brown   /* West */
311c4762a1bSJed Brown   h[0] = StaggeredMidpoint2D(p[i][j].h,p[i-1][j].h,p[i-1][j-1].h,p[i][j-1].h);
312c4762a1bSJed Brown   h[1] = StaggeredMidpoint2D(p[i][j].h,p[i-1][j].h,p[i-1][j+1].h,p[i][j+1].h);
313c4762a1bSJed Brown   /* East */
314c4762a1bSJed Brown   h[2] = StaggeredMidpoint2D(p[i][j].h,p[i+1][j].h,p[i+1][j+1].h,p[i][j+1].h);
315c4762a1bSJed Brown   h[3] = StaggeredMidpoint2D(p[i][j].h,p[i+1][j].h,p[i+1][j-1].h,p[i][j-1].h);
316c4762a1bSJed Brown   /* South */
317c4762a1bSJed Brown   h[4] = StaggeredMidpoint2D(p[i][j].h,p[i][j-1].h,p[i+1][j-1].h,p[i+1][j].h);
318c4762a1bSJed Brown   h[5] = StaggeredMidpoint2D(p[i][j].h,p[i][j-1].h,p[i-1][j-1].h,p[i-1][j].h);
319c4762a1bSJed Brown   /* North */
320c4762a1bSJed Brown   h[6] = StaggeredMidpoint2D(p[i][j].h,p[i][j+1].h,p[i-1][j+1].h,p[i-1][j].h);
321c4762a1bSJed Brown   h[7] = StaggeredMidpoint2D(p[i][j].h,p[i][j+1].h,p[i+1][j+1].h,p[i+1][j].h);
322c4762a1bSJed Brown }
323c4762a1bSJed Brown 
324c4762a1bSJed Brown /* Tests A and C are from the ISMIP-HOM paper (Pattyn et al. 2008) */
325c4762a1bSJed Brown static void THIInitialize_HOM_A(THI thi,PetscReal x,PetscReal y,PrmNode *p)
326c4762a1bSJed Brown {
327c4762a1bSJed Brown   Units units = thi->units;
328c4762a1bSJed Brown   PetscReal s = -x*PetscSinReal(thi->alpha);
329c4762a1bSJed Brown   p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx) * PetscSinReal(y*2*PETSC_PI/thi->Ly);
330c4762a1bSJed Brown   p->h = s - p->b;
331c4762a1bSJed 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  */
332c4762a1bSJed Brown }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown static void THIInitialize_HOM_C(THI thi,PetscReal x,PetscReal y,PrmNode *p)
335c4762a1bSJed Brown {
336c4762a1bSJed Brown   Units units = thi->units;
337c4762a1bSJed Brown   PetscReal s = -x*PetscSinReal(thi->alpha);
338c4762a1bSJed Brown   p->b = s - 1000*units->meter;
339c4762a1bSJed Brown   p->h = s - p->b;
340c4762a1bSJed Brown   /* tau_b = beta2 v   is a stress (Pa).
341c4762a1bSJed 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. */
342c4762a1bSJed 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;
343c4762a1bSJed Brown }
344c4762a1bSJed Brown 
345c4762a1bSJed Brown /* These are just toys */
346c4762a1bSJed Brown 
347c4762a1bSJed Brown /* From Fred Herman */
348c4762a1bSJed Brown static void THIInitialize_HOM_F(THI thi,PetscReal x,PetscReal y,PrmNode *p)
349c4762a1bSJed Brown {
350c4762a1bSJed Brown   Units units = thi->units;
351c4762a1bSJed Brown   PetscReal s = -x*PetscSinReal(thi->alpha);
352c4762a1bSJed Brown   p->b = s - 1000*units->meter + 100*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx);/* * sin(y*2*PETSC_PI/thi->Ly); */
353c4762a1bSJed Brown   p->h = s - p->b;
354c4762a1bSJed Brown   p->h = (1-(atan((x-thi->Lx/2)/1.)+PETSC_PI/2.)/PETSC_PI)*500*units->meter+1*units->meter;
355c4762a1bSJed Brown   s = PetscRealPart(p->b + p->h);
356c4762a1bSJed Brown   p->beta2 = -1e-10;
357c4762a1bSJed Brown   /*  p->beta2 = 1000 * units->Pascal * units->year / units->meter; */
358c4762a1bSJed Brown }
359c4762a1bSJed Brown 
360c4762a1bSJed Brown /* Same bed as test A, free slip everywhere except for a discontinuous jump to a circular sticky region in the middle. */
361c4762a1bSJed Brown static void THIInitialize_HOM_X(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
362c4762a1bSJed Brown {
363c4762a1bSJed Brown   Units units = thi->units;
364c4762a1bSJed Brown   PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
365c4762a1bSJed Brown   PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
366c4762a1bSJed Brown   p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
367c4762a1bSJed Brown   p->h = s - p->b;
368c4762a1bSJed Brown   p->beta2 = 1000 * (r < 1 ? 2 : 0) * units->Pascal * units->year / units->meter / thi->rhog;
369c4762a1bSJed Brown }
370c4762a1bSJed Brown 
371c4762a1bSJed Brown /* Like Z, but with 200 meter cliffs */
372c4762a1bSJed Brown static void THIInitialize_HOM_Y(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
373c4762a1bSJed Brown {
374c4762a1bSJed Brown   Units units = thi->units;
375c4762a1bSJed Brown   PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
376c4762a1bSJed Brown   PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
377c4762a1bSJed Brown   p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
378c4762a1bSJed Brown   if (PetscRealPart(p->b) > -700*units->meter) p->b += 200*units->meter;
379c4762a1bSJed Brown   p->h = s - p->b;
380c4762a1bSJed 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;
381c4762a1bSJed Brown }
382c4762a1bSJed Brown 
383c4762a1bSJed Brown /* Same bed as A, smoothly varying slipperiness, similar to MATLAB's "sombrero" (uncorrelated with bathymetry) */
384c4762a1bSJed Brown static void THIInitialize_HOM_Z(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
385c4762a1bSJed Brown {
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 
394c4762a1bSJed Brown static void THIFriction(THI thi,PetscReal rbeta2,PetscReal gam,PetscReal *beta2,PetscReal *dbeta2)
395c4762a1bSJed Brown {
396c4762a1bSJed Brown   if (thi->friction.irefgam == 0) {
397c4762a1bSJed Brown     Units units = thi->units;
398c4762a1bSJed Brown     thi->friction.irefgam = 1./(0.5*PetscSqr(100 * units->meter / units->year));
399c4762a1bSJed Brown     thi->friction.eps2 = 0.5*PetscSqr(1.e-4 / thi->friction.irefgam);
400c4762a1bSJed Brown   }
401c4762a1bSJed Brown   if (thi->friction.exponent == 0) {
402c4762a1bSJed Brown     *beta2  = rbeta2;
403c4762a1bSJed Brown     *dbeta2 = 0;
404c4762a1bSJed Brown   } else {
405c4762a1bSJed Brown     *beta2  = rbeta2 * PetscPowReal(thi->friction.eps2 + gam*thi->friction.irefgam,thi->friction.exponent);
406c4762a1bSJed Brown     *dbeta2 = thi->friction.exponent * *beta2 / (thi->friction.eps2 + gam*thi->friction.irefgam) * thi->friction.irefgam;
407c4762a1bSJed Brown   }
408c4762a1bSJed Brown }
409c4762a1bSJed Brown 
410c4762a1bSJed Brown static void THIViscosity(THI thi,PetscReal gam,PetscReal *eta,PetscReal *deta)
411c4762a1bSJed Brown {
412c4762a1bSJed Brown   PetscReal Bd2,eps,exponent;
413c4762a1bSJed Brown   if (thi->viscosity.Bd2 == 0) {
414c4762a1bSJed Brown     Units units = thi->units;
415c4762a1bSJed Brown     const PetscReal
416c4762a1bSJed Brown       n = thi->viscosity.glen_n,                        /* Glen exponent */
417c4762a1bSJed Brown       p = 1. + 1./n,                                    /* for Stokes */
418c4762a1bSJed Brown       A = 1.e-16 * PetscPowReal(units->Pascal,-n) / units->year, /* softness parameter (Pa^{-n}/s) */
419c4762a1bSJed Brown       B = PetscPowReal(A,-1./n);                                 /* hardness parameter */
420c4762a1bSJed Brown     thi->viscosity.Bd2      = B/2;
421c4762a1bSJed Brown     thi->viscosity.exponent = (p-2)/2;
422c4762a1bSJed Brown     thi->viscosity.eps      = 0.5*PetscSqr(1e-5 / units->year);
423c4762a1bSJed Brown   }
424c4762a1bSJed Brown   Bd2      = thi->viscosity.Bd2;
425c4762a1bSJed Brown   exponent = thi->viscosity.exponent;
426c4762a1bSJed Brown   eps      = thi->viscosity.eps;
427c4762a1bSJed Brown   *eta     = Bd2 * PetscPowReal(eps + gam,exponent);
428c4762a1bSJed Brown   *deta    = exponent * (*eta) / (eps + gam);
429c4762a1bSJed Brown }
430c4762a1bSJed Brown 
431c4762a1bSJed Brown static void THIErosion(THI thi,const Node *vel,PetscScalar *erate,Node *derate)
432c4762a1bSJed Brown {
433c4762a1bSJed Brown   const PetscScalar magref2 = 1.e-10 + (PetscSqr(vel->u) + PetscSqr(vel->v)) / PetscSqr(thi->erosion.refvel),
434c4762a1bSJed Brown                     rate    = -thi->erosion.rate*PetscPowScalar(magref2, 0.5*thi->erosion.exponent);
435c4762a1bSJed Brown   if (erate) *erate = rate;
436c4762a1bSJed Brown   if (derate) {
437c4762a1bSJed Brown     if (thi->erosion.exponent == 1) {
438c4762a1bSJed Brown       derate->u = 0;
439c4762a1bSJed Brown       derate->v = 0;
440c4762a1bSJed Brown     } else {
441c4762a1bSJed Brown       derate->u = 0.5*thi->erosion.exponent * rate / magref2 * 2. * vel->u / PetscSqr(thi->erosion.refvel);
442c4762a1bSJed Brown       derate->v = 0.5*thi->erosion.exponent * rate / magref2 * 2. * vel->v / PetscSqr(thi->erosion.refvel);
443c4762a1bSJed Brown     }
444c4762a1bSJed Brown   }
445c4762a1bSJed Brown }
446c4762a1bSJed Brown 
447c4762a1bSJed Brown static void RangeUpdate(PetscReal *min,PetscReal *max,PetscReal x)
448c4762a1bSJed Brown {
449c4762a1bSJed Brown   if (x < *min) *min = x;
450c4762a1bSJed Brown   if (x > *max) *max = x;
451c4762a1bSJed Brown }
452c4762a1bSJed Brown 
453c4762a1bSJed Brown static void PRangeClear(PRange *p)
454c4762a1bSJed Brown {
455c4762a1bSJed Brown   p->cmin = p->min = 1e100;
456c4762a1bSJed Brown   p->cmax = p->max = -1e100;
457c4762a1bSJed Brown }
458c4762a1bSJed Brown 
459c4762a1bSJed Brown static PetscErrorCode PRangeMinMax(PRange *p,PetscReal min,PetscReal max)
460c4762a1bSJed Brown {
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   PetscFunctionBeginUser;
463c4762a1bSJed Brown   p->cmin = min;
464c4762a1bSJed Brown   p->cmax = max;
465c4762a1bSJed Brown   if (min < p->min) p->min = min;
466c4762a1bSJed Brown   if (max > p->max) p->max = max;
467c4762a1bSJed Brown   PetscFunctionReturn(0);
468c4762a1bSJed Brown }
469c4762a1bSJed Brown 
470c4762a1bSJed Brown static PetscErrorCode THIDestroy(THI *thi)
471c4762a1bSJed Brown {
472c4762a1bSJed Brown   PetscErrorCode ierr;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   PetscFunctionBeginUser;
475c4762a1bSJed Brown   if (--((PetscObject)(*thi))->refct > 0) PetscFunctionReturn(0);
476c4762a1bSJed Brown   ierr = PetscFree((*thi)->units);CHKERRQ(ierr);
477c4762a1bSJed Brown   ierr = PetscFree((*thi)->mattype);CHKERRQ(ierr);
478c4762a1bSJed Brown   ierr = PetscFree((*thi)->monitor_basename);CHKERRQ(ierr);
479c4762a1bSJed Brown   ierr = PetscHeaderDestroy(thi);CHKERRQ(ierr);
480c4762a1bSJed Brown   PetscFunctionReturn(0);
481c4762a1bSJed Brown }
482c4762a1bSJed Brown 
483c4762a1bSJed Brown static PetscErrorCode THICreate(MPI_Comm comm,THI *inthi)
484c4762a1bSJed Brown {
485c4762a1bSJed Brown   static PetscBool registered = PETSC_FALSE;
486c4762a1bSJed Brown   THI              thi;
487c4762a1bSJed Brown   Units            units;
488c4762a1bSJed Brown   char             monitor_basename[PETSC_MAX_PATH_LEN] = "thi-";
489c4762a1bSJed Brown   PetscErrorCode   ierr;
490c4762a1bSJed Brown 
491c4762a1bSJed Brown   PetscFunctionBeginUser;
492c4762a1bSJed Brown   *inthi = 0;
493c4762a1bSJed Brown   if (!registered) {
494c4762a1bSJed Brown     ierr = PetscClassIdRegister("Toy Hydrostatic Ice",&THI_CLASSID);CHKERRQ(ierr);
495c4762a1bSJed Brown     registered = PETSC_TRUE;
496c4762a1bSJed Brown   }
497c4762a1bSJed Brown   ierr = PetscHeaderCreate(thi,THI_CLASSID,"THI","Toy Hydrostatic Ice","THI",comm,THIDestroy,0);CHKERRQ(ierr);
498c4762a1bSJed Brown 
499c4762a1bSJed Brown   ierr = PetscNew(&thi->units);CHKERRQ(ierr);
500c4762a1bSJed Brown 
501c4762a1bSJed Brown   units           = thi->units;
502c4762a1bSJed Brown   units->meter    = 1e-2;
503c4762a1bSJed Brown   units->second   = 1e-7;
504c4762a1bSJed Brown   units->kilogram = 1e-12;
505c4762a1bSJed Brown 
506c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"Scaled units options","");CHKERRQ(ierr);
507c4762a1bSJed Brown   {
508c4762a1bSJed Brown     ierr = PetscOptionsReal("-units_meter","1 meter in scaled length units","",units->meter,&units->meter,NULL);CHKERRQ(ierr);
509c4762a1bSJed Brown     ierr = PetscOptionsReal("-units_second","1 second in scaled time units","",units->second,&units->second,NULL);CHKERRQ(ierr);
510c4762a1bSJed Brown     ierr = PetscOptionsReal("-units_kilogram","1 kilogram in scaled mass units","",units->kilogram,&units->kilogram,NULL);CHKERRQ(ierr);
511c4762a1bSJed Brown   }
512c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
513c4762a1bSJed Brown   units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second));
514c4762a1bSJed Brown   units->year   = 31556926. * units->second, /* seconds per year */
515c4762a1bSJed Brown 
516c4762a1bSJed Brown   thi->Lx              = 10.e3;
517c4762a1bSJed Brown   thi->Ly              = 10.e3;
518c4762a1bSJed Brown   thi->Lz              = 1000;
519c4762a1bSJed Brown   thi->nlevels         = 1;
520c4762a1bSJed Brown   thi->dirichlet_scale = 1;
521c4762a1bSJed Brown   thi->verbose         = PETSC_FALSE;
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   thi->viscosity.glen_n = 3.;
524c4762a1bSJed Brown   thi->erosion.rate     = 1e-3; /* m/a */
525c4762a1bSJed Brown   thi->erosion.exponent = 1.;
526c4762a1bSJed Brown   thi->erosion.refvel   = 1.;   /* m/a */
527c4762a1bSJed Brown 
528c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"Toy Hydrostatic Ice options","");CHKERRQ(ierr);
529c4762a1bSJed Brown   {
530c4762a1bSJed Brown     QuadratureType quad       = QUAD_GAUSS;
531c4762a1bSJed Brown     char           homexp[]   = "A";
532c4762a1bSJed Brown     char           mtype[256] = MATSBAIJ;
533c4762a1bSJed Brown     PetscReal      L,m = 1.0;
534c4762a1bSJed Brown     PetscBool      flg;
535c4762a1bSJed Brown     L    = thi->Lx;
536c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_L","Domain size (m)","",L,&L,&flg);CHKERRQ(ierr);
537c4762a1bSJed Brown     if (flg) thi->Lx = thi->Ly = L;
538c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_Lx","X Domain size (m)","",thi->Lx,&thi->Lx,NULL);CHKERRQ(ierr);
539c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_Ly","Y Domain size (m)","",thi->Ly,&thi->Ly,NULL);CHKERRQ(ierr);
540c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_Lz","Z Domain size (m)","",thi->Lz,&thi->Lz,NULL);CHKERRQ(ierr);
541c4762a1bSJed Brown     ierr = PetscOptionsString("-thi_hom","ISMIP-HOM experiment (A or C)","",homexp,homexp,sizeof(homexp),NULL);CHKERRQ(ierr);
542c4762a1bSJed Brown     switch (homexp[0] = toupper(homexp[0])) {
543c4762a1bSJed Brown     case 'A':
544c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_A;
545c4762a1bSJed Brown       thi->no_slip    = PETSC_TRUE;
546c4762a1bSJed Brown       thi->alpha      = 0.5;
547c4762a1bSJed Brown       break;
548c4762a1bSJed Brown     case 'C':
549c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_C;
550c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
551c4762a1bSJed Brown       thi->alpha      = 0.1;
552c4762a1bSJed Brown       break;
553c4762a1bSJed Brown     case 'F':
554c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_F;
555c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
556c4762a1bSJed Brown       thi->alpha      = 0.5;
557c4762a1bSJed Brown       break;
558c4762a1bSJed Brown     case 'X':
559c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_X;
560c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
561c4762a1bSJed Brown       thi->alpha      = 0.3;
562c4762a1bSJed Brown       break;
563c4762a1bSJed Brown     case 'Y':
564c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_Y;
565c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
566c4762a1bSJed Brown       thi->alpha      = 0.5;
567c4762a1bSJed Brown       break;
568c4762a1bSJed Brown     case 'Z':
569c4762a1bSJed Brown       thi->initialize = THIInitialize_HOM_Z;
570c4762a1bSJed Brown       thi->no_slip    = PETSC_FALSE;
571c4762a1bSJed Brown       thi->alpha      = 0.5;
572c4762a1bSJed Brown       break;
573c4762a1bSJed Brown     default:
574c4762a1bSJed Brown       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"HOM experiment '%c' not implemented",homexp[0]);
575c4762a1bSJed Brown     }
576c4762a1bSJed Brown     ierr = PetscOptionsEnum("-thi_quadrature","Quadrature to use for 3D elements","",QuadratureTypes,(PetscEnum)quad,(PetscEnum*)&quad,NULL);CHKERRQ(ierr);
577c4762a1bSJed Brown     switch (quad) {
578c4762a1bSJed Brown     case QUAD_GAUSS:
579c4762a1bSJed Brown       HexQInterp = HexQInterp_Gauss;
580c4762a1bSJed Brown       HexQDeriv  = HexQDeriv_Gauss;
581c4762a1bSJed Brown       break;
582c4762a1bSJed Brown     case QUAD_LOBATTO:
583c4762a1bSJed Brown       HexQInterp = HexQInterp_Lobatto;
584c4762a1bSJed Brown       HexQDeriv  = HexQDeriv_Lobatto;
585c4762a1bSJed Brown       break;
586c4762a1bSJed Brown     }
587c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_alpha","Bed angle (degrees)","",thi->alpha,&thi->alpha,NULL);CHKERRQ(ierr);
588c4762a1bSJed Brown     ierr = 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);CHKERRQ(ierr);
589c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_friction_m","Friction exponent, 0=Coulomb, 1=Navier","",m,&m,NULL);CHKERRQ(ierr);
590c4762a1bSJed Brown     thi->friction.exponent = (m-1)/2;
591c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_erosion_rate","Rate of erosion relative to sliding velocity at reference velocity (m/a)",NULL,thi->erosion.rate,&thi->erosion.rate,NULL);CHKERRQ(ierr);
592c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_erosion_exponent","Power of sliding velocity appearing in erosion relation",NULL,thi->erosion.exponent,&thi->erosion.exponent,NULL);CHKERRQ(ierr);
593c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_erosion_refvel","Reference sliding velocity for erosion (m/a)",NULL,thi->erosion.refvel,&thi->erosion.refvel,NULL);CHKERRQ(ierr);
594c4762a1bSJed Brown     thi->erosion.rate   *= units->meter / units->year;
595c4762a1bSJed Brown     thi->erosion.refvel *= units->meter / units->year;
596c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_dirichlet_scale","Scale Dirichlet boundary conditions by this factor","",thi->dirichlet_scale,&thi->dirichlet_scale,NULL);CHKERRQ(ierr);
597c4762a1bSJed Brown     ierr = 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);CHKERRQ(ierr);
598c4762a1bSJed Brown     ierr = PetscOptionsReal("-thi_inertia","Coefficient of accelaration term in velocity system, physical is almost zero",NULL,thi->inertia,&thi->inertia,NULL);CHKERRQ(ierr);
599c4762a1bSJed Brown     ierr = PetscOptionsInt("-thi_nlevels","Number of levels of refinement","",thi->nlevels,&thi->nlevels,NULL);CHKERRQ(ierr);
600c4762a1bSJed Brown     ierr = PetscOptionsFList("-thi_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)mtype,sizeof(mtype),NULL);CHKERRQ(ierr);
601c4762a1bSJed Brown     ierr = PetscStrallocpy(mtype,&thi->mattype);CHKERRQ(ierr);
602c4762a1bSJed Brown     ierr = PetscOptionsBool("-thi_verbose","Enable verbose output (like matrix sizes and statistics)","",thi->verbose,&thi->verbose,NULL);CHKERRQ(ierr);
603c4762a1bSJed Brown     ierr = PetscOptionsString("-thi_monitor","Basename to write state files to",NULL,monitor_basename,monitor_basename,sizeof(monitor_basename),&flg);CHKERRQ(ierr);
604c4762a1bSJed Brown     if (flg) {
605c4762a1bSJed Brown       ierr = PetscStrallocpy(monitor_basename,&thi->monitor_basename);CHKERRQ(ierr);
606c4762a1bSJed Brown       thi->monitor_interval = 1;
607c4762a1bSJed Brown       ierr = PetscOptionsInt("-thi_monitor_interval","Frequency at which to write state files",NULL,thi->monitor_interval,&thi->monitor_interval,NULL);CHKERRQ(ierr);
608c4762a1bSJed Brown     }
609c4762a1bSJed Brown   }
610c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
611c4762a1bSJed Brown 
612c4762a1bSJed Brown   /* dimensionalize */
613c4762a1bSJed Brown   thi->Lx    *= units->meter;
614c4762a1bSJed Brown   thi->Ly    *= units->meter;
615c4762a1bSJed Brown   thi->Lz    *= units->meter;
616c4762a1bSJed Brown   thi->alpha *= PETSC_PI / 180;
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   PRangeClear(&thi->eta);
619c4762a1bSJed Brown   PRangeClear(&thi->beta2);
620c4762a1bSJed Brown 
621c4762a1bSJed Brown   {
622c4762a1bSJed Brown     PetscReal u       = 1000*units->meter/(3e7*units->second),
623c4762a1bSJed Brown               gradu   = u / (100*units->meter),eta,deta,
624c4762a1bSJed Brown               rho     = 910 * units->kilogram/PetscPowRealInt(units->meter,3),
625c4762a1bSJed Brown               grav    = 9.81 * units->meter/PetscSqr(units->second),
626c4762a1bSJed Brown               driving = rho * grav * PetscSinReal(thi->alpha) * 1000*units->meter;
627c4762a1bSJed Brown     THIViscosity(thi,0.5*gradu*gradu,&eta,&deta);
628c4762a1bSJed Brown     thi->rhog = rho * grav;
629c4762a1bSJed Brown     if (thi->verbose) {
630c4762a1bSJed Brown       ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Units: meter %8.2g  second %8.2g  kg %8.2g  Pa %8.2g\n",units->meter,units->second,units->kilogram,units->Pascal);CHKERRQ(ierr);
631c4762a1bSJed Brown       ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Domain (%6.2g,%6.2g,%6.2g), pressure %8.2g, driving stress %8.2g\n",thi->Lx,thi->Ly,thi->Lz,rho*grav*1e3*units->meter,driving);CHKERRQ(ierr);
632c4762a1bSJed Brown       ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Large velocity 1km/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",u,gradu,eta,2*eta*gradu,2*eta*gradu/driving);CHKERRQ(ierr);
633c4762a1bSJed Brown       THIViscosity(thi,0.5*PetscSqr(1e-3*gradu),&eta,&deta);
634c4762a1bSJed Brown       ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Small velocity 1m/a  %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",1e-3*u,1e-3*gradu,eta,2*eta*1e-3*gradu,2*eta*1e-3*gradu/driving);CHKERRQ(ierr);
635c4762a1bSJed Brown     }
636c4762a1bSJed Brown   }
637c4762a1bSJed Brown 
638c4762a1bSJed Brown   *inthi = thi;
639c4762a1bSJed Brown   PetscFunctionReturn(0);
640c4762a1bSJed Brown }
641c4762a1bSJed Brown 
642c4762a1bSJed Brown /* Our problem is periodic, but the domain has a mean slope of alpha so the bed does not line up between the upstream
643c4762a1bSJed Brown  * and downstream ends of the domain.  This function fixes the ghost values so that the domain appears truly periodic in
644c4762a1bSJed Brown  * the horizontal. */
645c4762a1bSJed Brown static PetscErrorCode THIFixGhosts(THI thi,DM da3,DM da2,Vec X3,Vec X2)
646c4762a1bSJed Brown {
647c4762a1bSJed Brown   PetscErrorCode ierr;
648c4762a1bSJed Brown   DMDALocalInfo  info;
649c4762a1bSJed Brown   PrmNode        **x2;
650c4762a1bSJed Brown   PetscInt       i,j;
651c4762a1bSJed Brown 
652c4762a1bSJed Brown   PetscFunctionBeginUser;
653c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da3,&info);CHKERRQ(ierr);
654c4762a1bSJed Brown   /* ierr = VecView(X2,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
655c4762a1bSJed Brown   ierr = DMDAVecGetArray(da2,X2,&x2);CHKERRQ(ierr);
656c4762a1bSJed Brown   for (i=info.gzs; i<info.gzs+info.gzm; i++) {
657c4762a1bSJed Brown     if (i > -1 && i < info.mz) continue;
658c4762a1bSJed Brown     for (j=info.gys; j<info.gys+info.gym; j++) {
659c4762a1bSJed Brown       x2[i][j].b += PetscSinReal(thi->alpha) * thi->Lx * (i<0 ? 1.0 : -1.0);
660c4762a1bSJed Brown     }
661c4762a1bSJed Brown   }
662c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da2,X2,&x2);CHKERRQ(ierr);
663c4762a1bSJed Brown   /* ierr = VecView(X2,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
664c4762a1bSJed Brown   PetscFunctionReturn(0);
665c4762a1bSJed Brown }
666c4762a1bSJed Brown 
667c4762a1bSJed Brown static PetscErrorCode THIInitializePrm(THI thi,DM da2prm,PrmNode **p)
668c4762a1bSJed Brown {
669c4762a1bSJed Brown   PetscInt       i,j,xs,xm,ys,ym,mx,my;
670c4762a1bSJed Brown   PetscErrorCode ierr;
671c4762a1bSJed Brown 
672c4762a1bSJed Brown   PetscFunctionBeginUser;
673c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(da2prm,&ys,&xs,0,&ym,&xm,0);CHKERRQ(ierr);
674c4762a1bSJed Brown   ierr = DMDAGetInfo(da2prm,0, &my,&mx,0, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
675c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
676c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
677c4762a1bSJed Brown       PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my;
678c4762a1bSJed Brown       thi->initialize(thi,xx,yy,&p[i][j]);
679c4762a1bSJed Brown     }
680c4762a1bSJed Brown   }
681c4762a1bSJed Brown   PetscFunctionReturn(0);
682c4762a1bSJed Brown }
683c4762a1bSJed Brown 
684c4762a1bSJed Brown static PetscErrorCode THIInitial(THI thi,DM pack,Vec X)
685c4762a1bSJed Brown {
686c4762a1bSJed Brown   DM             da3,da2;
687c4762a1bSJed Brown   PetscInt       i,j,k,xs,xm,ys,ym,zs,zm,mx,my;
688c4762a1bSJed Brown   PetscReal      hx,hy;
689c4762a1bSJed Brown   PrmNode        **prm;
690c4762a1bSJed Brown   Node           ***x;
691c4762a1bSJed Brown   Vec            X3g,X2g,X2;
692c4762a1bSJed Brown   PetscErrorCode ierr;
693c4762a1bSJed Brown 
694c4762a1bSJed Brown   PetscFunctionBeginUser;
695c4762a1bSJed Brown   ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr);
696c4762a1bSJed Brown   ierr = DMCompositeGetAccess(pack,X,&X3g,&X2g);CHKERRQ(ierr);
697c4762a1bSJed Brown   ierr = DMGetLocalVector(da2,&X2);CHKERRQ(ierr);
698c4762a1bSJed Brown 
699c4762a1bSJed Brown   ierr = DMDAGetInfo(da3,0, 0,&my,&mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
700c4762a1bSJed Brown   ierr = DMDAGetCorners(da3,&zs,&ys,&xs,&zm,&ym,&xm);CHKERRQ(ierr);
701c4762a1bSJed Brown   ierr = DMDAVecGetArray(da3,X3g,&x);CHKERRQ(ierr);
702c4762a1bSJed Brown   ierr = DMDAVecGetArray(da2,X2,&prm);CHKERRQ(ierr);
703c4762a1bSJed Brown 
704c4762a1bSJed Brown   ierr = THIInitializePrm(thi,da2,prm);CHKERRQ(ierr);
705c4762a1bSJed Brown 
706c4762a1bSJed Brown   hx = thi->Lx / mx;
707c4762a1bSJed Brown   hy = thi->Ly / my;
708c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
709c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
710c4762a1bSJed Brown       for (k=zs; k<zs+zm; k++) {
711c4762a1bSJed Brown         const PetscScalar zm1      = zm-1,
712c4762a1bSJed Brown                           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),
713c4762a1bSJed Brown                           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);
714c4762a1bSJed Brown         x[i][j][k].u = 0. * drivingx * prm[i][j].h*(PetscScalar)k/zm1;
715c4762a1bSJed Brown         x[i][j][k].v = 0. * drivingy * prm[i][j].h*(PetscScalar)k/zm1;
716c4762a1bSJed Brown       }
717c4762a1bSJed Brown     }
718c4762a1bSJed Brown   }
719c4762a1bSJed Brown 
720c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da3,X3g,&x);CHKERRQ(ierr);
721c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da2,X2,&prm);CHKERRQ(ierr);
722c4762a1bSJed Brown 
723c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da2,X2,INSERT_VALUES,X2g);CHKERRQ(ierr);
724c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd  (da2,X2,INSERT_VALUES,X2g);CHKERRQ(ierr);
725c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da2,&X2);CHKERRQ(ierr);
726c4762a1bSJed Brown 
727c4762a1bSJed Brown   ierr = DMCompositeRestoreAccess(pack,X,&X3g,&X2g);CHKERRQ(ierr);
728c4762a1bSJed Brown   PetscFunctionReturn(0);
729c4762a1bSJed Brown }
730c4762a1bSJed Brown 
731c4762a1bSJed Brown 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)
732c4762a1bSJed Brown {
733c4762a1bSJed Brown   PetscInt    l,ll;
734c4762a1bSJed Brown   PetscScalar gam;
735c4762a1bSJed Brown 
736c4762a1bSJed Brown   du[0] = du[1] = du[2] = 0;
737c4762a1bSJed Brown   dv[0] = dv[1] = dv[2] = 0;
738c4762a1bSJed Brown   *u    = 0;
739c4762a1bSJed Brown   *v    = 0;
740c4762a1bSJed Brown   for (l=0; l<8; l++) {
741c4762a1bSJed Brown     *u += phi[l] * n[l].u;
742c4762a1bSJed Brown     *v += phi[l] * n[l].v;
743c4762a1bSJed Brown     for (ll=0; ll<3; ll++) {
744c4762a1bSJed Brown       du[ll] += dphi[l][ll] * n[l].u;
745c4762a1bSJed Brown       dv[ll] += dphi[l][ll] * n[l].v;
746c4762a1bSJed Brown     }
747c4762a1bSJed Brown   }
748c4762a1bSJed 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]);
749c4762a1bSJed Brown   THIViscosity(thi,PetscRealPart(gam),eta,deta);
750c4762a1bSJed Brown }
751c4762a1bSJed Brown 
752c4762a1bSJed Brown static PetscErrorCode THIFunctionLocal_3D(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,const Node ***xdot,Node ***f,THI thi)
753c4762a1bSJed Brown {
754c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,zm,i,j,k,q,l;
755c4762a1bSJed Brown   PetscReal      hx,hy,etamin,etamax,beta2min,beta2max;
756c4762a1bSJed Brown   PetscErrorCode ierr;
757c4762a1bSJed Brown 
758c4762a1bSJed Brown   PetscFunctionBeginUser;
759c4762a1bSJed Brown   xs = info->zs;
760c4762a1bSJed Brown   ys = info->ys;
761c4762a1bSJed Brown   xm = info->zm;
762c4762a1bSJed Brown   ym = info->ym;
763c4762a1bSJed Brown   zm = info->xm;
764c4762a1bSJed Brown   hx = thi->Lx / info->mz;
765c4762a1bSJed Brown   hy = thi->Ly / info->my;
766c4762a1bSJed Brown 
767c4762a1bSJed Brown   etamin   = 1e100;
768c4762a1bSJed Brown   etamax   = 0;
769c4762a1bSJed Brown   beta2min = 1e100;
770c4762a1bSJed Brown   beta2max = 0;
771c4762a1bSJed Brown 
772c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
773c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
774c4762a1bSJed Brown       PrmNode pn[4],dpn[4][2];
775c4762a1bSJed Brown       QuadExtract(prm,i,j,pn);
776c4762a1bSJed Brown       ierr = QuadComputeGrad4(QuadQDeriv,hx,hy,pn,dpn);CHKERRQ(ierr);
777c4762a1bSJed Brown       for (k=0; k<zm-1; k++) {
778c4762a1bSJed Brown         PetscInt  ls = 0;
779c4762a1bSJed Brown         Node      n[8],ndot[8],*fn[8];
780c4762a1bSJed Brown         PetscReal zn[8],etabase = 0;
781*2f613bf5SBarry Smith 
782c4762a1bSJed Brown         PrmHexGetZ(pn,k,zm,zn);
783c4762a1bSJed Brown         HexExtract(x,i,j,k,n);
784*2f613bf5SBarry Smith         HexExtract(xdot,i,j,k,ndot);
785c4762a1bSJed Brown         HexExtractRef(f,i,j,k,fn);
786c4762a1bSJed Brown         if (thi->no_slip && k == 0) {
787c4762a1bSJed Brown           for (l=0; l<4; l++) n[l].u = n[l].v = 0;
788c4762a1bSJed Brown           /* The first 4 basis functions lie on the bottom layer, so their contribution is exactly 0, hence we can skip them */
789c4762a1bSJed Brown           ls = 4;
790c4762a1bSJed Brown         }
791c4762a1bSJed Brown         for (q=0; q<8; q++) {
792c4762a1bSJed Brown           PetscReal   dz[3],phi[8],dphi[8][3],jw,eta,deta;
793c4762a1bSJed Brown           PetscScalar du[3],dv[3],u,v,udot=0,vdot=0;
794c4762a1bSJed Brown           for (l=ls; l<8; l++) {
795c4762a1bSJed Brown             udot += HexQInterp[q][l]*ndot[l].u;
796c4762a1bSJed Brown             vdot += HexQInterp[q][l]*ndot[l].v;
797c4762a1bSJed Brown           }
798c4762a1bSJed Brown           HexGrad(HexQDeriv[q],zn,dz);
799c4762a1bSJed Brown           HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
800c4762a1bSJed Brown           PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
801c4762a1bSJed Brown           jw /= thi->rhog;      /* scales residuals to be O(1) */
802c4762a1bSJed Brown           if (q == 0) etabase = eta;
803c4762a1bSJed Brown           RangeUpdate(&etamin,&etamax,eta);
804c4762a1bSJed Brown           for (l=ls; l<8; l++) { /* test functions */
805c4762a1bSJed 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};
806c4762a1bSJed Brown             const PetscReal   pp    = phi[l],*dp = dphi[l];
807c4762a1bSJed 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];
808c4762a1bSJed 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];
809c4762a1bSJed Brown             fn[l]->u += pp*jw*udot*thi->inertia*pp;
810c4762a1bSJed Brown             fn[l]->v += pp*jw*vdot*thi->inertia*pp;
811c4762a1bSJed Brown           }
812c4762a1bSJed Brown         }
813c4762a1bSJed Brown         if (k == 0) { /* we are on a bottom face */
814c4762a1bSJed Brown           if (thi->no_slip) {
815c4762a1bSJed Brown             /* Note: Non-Galerkin coarse grid operators are very sensitive to the scaling of Dirichlet boundary
816c4762a1bSJed Brown             * conditions.  After shenanigans above, etabase contains the effective viscosity at the closest quadrature
817c4762a1bSJed Brown             * point to the bed.  We want the diagonal entry in the Dirichlet condition to have similar magnitude to the
818c4762a1bSJed Brown             * diagonal entry corresponding to the adjacent node.  The fundamental scaling of the viscous part is in
819c4762a1bSJed Brown             * diagu, diagv below.  This scaling is easy to recognize by considering the finite difference operator after
820c4762a1bSJed Brown             * scaling by element size.  The no-slip Dirichlet condition is scaled by this factor, and also in the
821c4762a1bSJed Brown             * assembled matrix (see the similar block in THIJacobianLocal).
822c4762a1bSJed Brown             *
823c4762a1bSJed Brown             * Note that the residual at this Dirichlet node is linear in the state at this node, but also depends
824c4762a1bSJed Brown             * (nonlinearly in general) on the neighboring interior nodes through the local viscosity.  This will make
825c4762a1bSJed Brown             * a matrix-free Jacobian have extra entries in the corresponding row.  We assemble only the diagonal part,
826c4762a1bSJed Brown             * so the solution will exactly satisfy the boundary condition after the first linear iteration.
827c4762a1bSJed Brown             */
828c4762a1bSJed Brown             const PetscReal   hz    = PetscRealPart(pn[0].h)/(zm-1.);
829c4762a1bSJed 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);
830c4762a1bSJed Brown             fn[0]->u = thi->dirichlet_scale*diagu*x[i][j][k].u;
831c4762a1bSJed Brown             fn[0]->v = thi->dirichlet_scale*diagv*x[i][j][k].v;
832c4762a1bSJed Brown           } else {              /* Integrate over bottom face to apply boundary condition */
833c4762a1bSJed 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) */
834c4762a1bSJed Brown               const PetscReal jw = 0.25*hx*hy,*phi = QuadQInterp[q];
835c4762a1bSJed Brown               PetscScalar     u  =0,v=0,rbeta2=0;
836c4762a1bSJed Brown               PetscReal       beta2,dbeta2;
837c4762a1bSJed Brown               for (l=0; l<4; l++) {
838c4762a1bSJed Brown                 u     += phi[l]*n[l].u;
839c4762a1bSJed Brown                 v     += phi[l]*n[l].v;
840c4762a1bSJed Brown                 rbeta2 += phi[l]*pn[l].beta2;
841c4762a1bSJed Brown               }
842c4762a1bSJed Brown               THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
843c4762a1bSJed Brown               RangeUpdate(&beta2min,&beta2max,beta2);
844c4762a1bSJed Brown               for (l=0; l<4; l++) {
845c4762a1bSJed Brown                 const PetscReal pp = phi[l];
846c4762a1bSJed Brown                 fn[ls+l]->u += pp*jw*beta2*u;
847c4762a1bSJed Brown                 fn[ls+l]->v += pp*jw*beta2*v;
848c4762a1bSJed Brown               }
849c4762a1bSJed Brown             }
850c4762a1bSJed Brown           }
851c4762a1bSJed Brown         }
852c4762a1bSJed Brown       }
853c4762a1bSJed Brown     }
854c4762a1bSJed Brown   }
855c4762a1bSJed Brown 
856c4762a1bSJed Brown   ierr = PRangeMinMax(&thi->eta,etamin,etamax);CHKERRQ(ierr);
857c4762a1bSJed Brown   ierr = PRangeMinMax(&thi->beta2,beta2min,beta2max);CHKERRQ(ierr);
858c4762a1bSJed Brown   PetscFunctionReturn(0);
859c4762a1bSJed Brown }
860c4762a1bSJed Brown 
861c4762a1bSJed Brown static PetscErrorCode THIFunctionLocal_2D(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,const PrmNode **prmdot,PrmNode **f,THI thi)
862c4762a1bSJed Brown {
863c4762a1bSJed Brown   PetscInt xs,ys,xm,ym,zm,i,j,k;
864c4762a1bSJed Brown 
865c4762a1bSJed Brown   PetscFunctionBeginUser;
866c4762a1bSJed Brown   xs = info->zs;
867c4762a1bSJed Brown   ys = info->ys;
868c4762a1bSJed Brown   xm = info->zm;
869c4762a1bSJed Brown   ym = info->ym;
870c4762a1bSJed Brown   zm = info->xm;
871c4762a1bSJed Brown 
872c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
873c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
874c4762a1bSJed Brown       PetscScalar div = 0,erate,h[8];
875c4762a1bSJed Brown       PrmNodeGetFaceMeasure(prm,i,j,h);
876c4762a1bSJed Brown       for (k=0; k<zm; k++) {
877c4762a1bSJed Brown         PetscScalar weight = (k==0 || k == zm-1) ? 0.5/(zm-1) : 1.0/(zm-1);
878c4762a1bSJed Brown         if (0) {                /* centered flux */
879c4762a1bSJed Brown           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)
880c4762a1bSJed Brown                   - 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)
881c4762a1bSJed Brown                   + 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)
882c4762a1bSJed Brown                   + 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)
883c4762a1bSJed Brown                   - 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)
884c4762a1bSJed Brown                   - 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)
885c4762a1bSJed Brown                   + 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)
886c4762a1bSJed Brown                   + 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));
887c4762a1bSJed Brown         } else {                /* Upwind flux */
888c4762a1bSJed Brown           div += weight*(-UpwindFluxXW(x,prm,h,i,j,k, 1)
889c4762a1bSJed Brown                          -UpwindFluxXW(x,prm,h,i,j,k,-1)
890c4762a1bSJed Brown                          +UpwindFluxXE(x,prm,h,i,j,k, 1)
891c4762a1bSJed Brown                          +UpwindFluxXE(x,prm,h,i,j,k,-1)
892c4762a1bSJed Brown                          -UpwindFluxYS(x,prm,h,i,j,k, 1)
893c4762a1bSJed Brown                          -UpwindFluxYS(x,prm,h,i,j,k,-1)
894c4762a1bSJed Brown                          +UpwindFluxYN(x,prm,h,i,j,k, 1)
895c4762a1bSJed Brown                          +UpwindFluxYN(x,prm,h,i,j,k,-1));
896c4762a1bSJed Brown         }
897c4762a1bSJed Brown       }
898c4762a1bSJed Brown       /* printf("div[%d][%d] %g\n",i,j,div); */
899c4762a1bSJed Brown       THIErosion(thi,&x[i][j][0],&erate,NULL);
900c4762a1bSJed Brown       f[i][j].b     = prmdot[i][j].b - erate;
901c4762a1bSJed Brown       f[i][j].h     = prmdot[i][j].h + div;
902c4762a1bSJed Brown       f[i][j].beta2 = prmdot[i][j].beta2;
903c4762a1bSJed Brown     }
904c4762a1bSJed Brown   }
905c4762a1bSJed Brown   PetscFunctionReturn(0);
906c4762a1bSJed Brown }
907c4762a1bSJed Brown 
908c4762a1bSJed Brown static PetscErrorCode THIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
909c4762a1bSJed Brown {
910c4762a1bSJed Brown   PetscErrorCode ierr;
911c4762a1bSJed Brown   THI            thi = (THI)ctx;
912c4762a1bSJed Brown   DM             pack,da3,da2;
913c4762a1bSJed Brown   Vec            X3,X2,Xdot3,Xdot2,F3,F2,F3g,F2g;
914c4762a1bSJed Brown   const Node     ***x3,***xdot3;
915c4762a1bSJed Brown   const PrmNode  **x2,**xdot2;
916c4762a1bSJed Brown   Node           ***f3;
917c4762a1bSJed Brown   PrmNode        **f2;
918c4762a1bSJed Brown   DMDALocalInfo  info3;
919c4762a1bSJed Brown 
920c4762a1bSJed Brown   PetscFunctionBeginUser;
921c4762a1bSJed Brown   ierr = TSGetDM(ts,&pack);CHKERRQ(ierr);
922c4762a1bSJed Brown   ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr);
923c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da3,&info3);CHKERRQ(ierr);
924c4762a1bSJed Brown   ierr = DMCompositeGetLocalVectors(pack,&X3,&X2);CHKERRQ(ierr);
925c4762a1bSJed Brown   ierr = DMCompositeGetLocalVectors(pack,&Xdot3,&Xdot2);CHKERRQ(ierr);
926c4762a1bSJed Brown   ierr = DMCompositeScatter(pack,X,X3,X2);CHKERRQ(ierr);
927c4762a1bSJed Brown   ierr = THIFixGhosts(thi,da3,da2,X3,X2);CHKERRQ(ierr);
928c4762a1bSJed Brown   ierr = DMCompositeScatter(pack,Xdot,Xdot3,Xdot2);CHKERRQ(ierr);
929c4762a1bSJed Brown 
930c4762a1bSJed Brown   ierr = DMGetLocalVector(da3,&F3);CHKERRQ(ierr);
931c4762a1bSJed Brown   ierr = DMGetLocalVector(da2,&F2);CHKERRQ(ierr);
932c4762a1bSJed Brown   ierr = VecZeroEntries(F3);CHKERRQ(ierr);
933c4762a1bSJed Brown 
934c4762a1bSJed Brown   ierr = DMDAVecGetArray(da3,X3,&x3);CHKERRQ(ierr);
935c4762a1bSJed Brown   ierr = DMDAVecGetArray(da2,X2,&x2);CHKERRQ(ierr);
936c4762a1bSJed Brown   ierr = DMDAVecGetArray(da3,Xdot3,&xdot3);CHKERRQ(ierr);
937c4762a1bSJed Brown   ierr = DMDAVecGetArray(da2,Xdot2,&xdot2);CHKERRQ(ierr);
938c4762a1bSJed Brown   ierr = DMDAVecGetArray(da3,F3,&f3);CHKERRQ(ierr);
939c4762a1bSJed Brown   ierr = DMDAVecGetArray(da2,F2,&f2);CHKERRQ(ierr);
940c4762a1bSJed Brown 
941c4762a1bSJed Brown   ierr = THIFunctionLocal_3D(&info3,x3,x2,xdot3,f3,thi);CHKERRQ(ierr);
942c4762a1bSJed Brown   ierr = THIFunctionLocal_2D(&info3,x3,x2,xdot2,f2,thi);CHKERRQ(ierr);
943c4762a1bSJed Brown 
944c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da3,X3,&x3);CHKERRQ(ierr);
945c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da2,X2,&x2);CHKERRQ(ierr);
946c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da3,Xdot3,&xdot3);CHKERRQ(ierr);
947c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da2,Xdot2,&xdot2);CHKERRQ(ierr);
948c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da3,F3,&f3);CHKERRQ(ierr);
949c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da2,F2,&f2);CHKERRQ(ierr);
950c4762a1bSJed Brown 
951c4762a1bSJed Brown   ierr = DMCompositeRestoreLocalVectors(pack,&X3,&X2);CHKERRQ(ierr);
952c4762a1bSJed Brown   ierr = DMCompositeRestoreLocalVectors(pack,&Xdot3,&Xdot2);CHKERRQ(ierr);
953c4762a1bSJed Brown 
954c4762a1bSJed Brown   ierr = VecZeroEntries(F);CHKERRQ(ierr);
955c4762a1bSJed Brown   ierr = DMCompositeGetAccess(pack,F,&F3g,&F2g);CHKERRQ(ierr);
956c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da3,F3,ADD_VALUES,F3g);CHKERRQ(ierr);
957c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd  (da3,F3,ADD_VALUES,F3g);CHKERRQ(ierr);
958c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da2,F2,INSERT_VALUES,F2g);CHKERRQ(ierr);
959c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd  (da2,F2,INSERT_VALUES,F2g);CHKERRQ(ierr);
960c4762a1bSJed Brown 
961c4762a1bSJed Brown   if (thi->verbose) {
962c4762a1bSJed Brown     PetscViewer viewer;
963c4762a1bSJed Brown     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)thi),&viewer);CHKERRQ(ierr);
964c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"3D_Velocity residual (bs=2):\n");CHKERRQ(ierr);
965c4762a1bSJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
966c4762a1bSJed Brown     ierr = VecView(F3,viewer);CHKERRQ(ierr);
967c4762a1bSJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
968c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"2D_Fields residual (bs=3):\n");CHKERRQ(ierr);
969c4762a1bSJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
970c4762a1bSJed Brown     ierr = VecView(F2,viewer);CHKERRQ(ierr);
971c4762a1bSJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
972c4762a1bSJed Brown   }
973c4762a1bSJed Brown 
974c4762a1bSJed Brown   ierr = DMCompositeRestoreAccess(pack,F,&F3g,&F2g);CHKERRQ(ierr);
975c4762a1bSJed Brown 
976c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da3,&F3);CHKERRQ(ierr);
977c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da2,&F2);CHKERRQ(ierr);
978c4762a1bSJed Brown   PetscFunctionReturn(0);
979c4762a1bSJed Brown }
980c4762a1bSJed Brown 
981c4762a1bSJed Brown static PetscErrorCode THIMatrixStatistics(THI thi,Mat B,PetscViewer viewer)
982c4762a1bSJed Brown {
983c4762a1bSJed Brown   PetscErrorCode ierr;
984c4762a1bSJed Brown   PetscReal      nrm;
985c4762a1bSJed Brown   PetscInt       m;
986c4762a1bSJed Brown   PetscMPIInt    rank;
987c4762a1bSJed Brown 
988c4762a1bSJed Brown   PetscFunctionBeginUser;
989c4762a1bSJed Brown   ierr = MatNorm(B,NORM_FROBENIUS,&nrm);CHKERRQ(ierr);
990c4762a1bSJed Brown   ierr = MatGetSize(B,&m,0);CHKERRQ(ierr);
991ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)B),&rank);CHKERRMPI(ierr);
992dd400576SPatrick Sanan   if (rank == 0) {
993c4762a1bSJed Brown     PetscScalar val0,val2;
994c4762a1bSJed Brown     ierr = MatGetValue(B,0,0,&val0);CHKERRQ(ierr);
995c4762a1bSJed Brown     ierr = MatGetValue(B,2,2,&val2);CHKERRQ(ierr);
996c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"Matrix dim %8d  norm %8.2e, (0,0) %8.2e  (2,2) %8.2e, eta [%8.2e,%8.2e] beta2 [%8.2e,%8.2e]\n",m,nrm,PetscRealPart(val0),PetscRealPart(val2),thi->eta.cmin,thi->eta.cmax,thi->beta2.cmin,thi->beta2.cmax);CHKERRQ(ierr);
997c4762a1bSJed Brown   }
998c4762a1bSJed Brown   PetscFunctionReturn(0);
999c4762a1bSJed Brown }
1000c4762a1bSJed Brown 
1001c4762a1bSJed Brown static PetscErrorCode THISurfaceStatistics(DM pack,Vec X,PetscReal *min,PetscReal *max,PetscReal *mean)
1002c4762a1bSJed Brown {
1003c4762a1bSJed Brown   PetscErrorCode ierr;
1004c4762a1bSJed Brown   DM             da3,da2;
1005c4762a1bSJed Brown   Vec            X3,X2;
1006c4762a1bSJed Brown   Node           ***x;
1007c4762a1bSJed Brown   PetscInt       i,j,xs,ys,zs,xm,ym,zm,mx,my,mz;
1008c4762a1bSJed Brown   PetscReal      umin = 1e100,umax=-1e100;
1009c4762a1bSJed Brown   PetscScalar    usum =0.0,gusum;
1010c4762a1bSJed Brown 
1011c4762a1bSJed Brown   PetscFunctionBeginUser;
1012c4762a1bSJed Brown   ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr);
1013c4762a1bSJed Brown   ierr = DMCompositeGetAccess(pack,X,&X3,&X2);CHKERRQ(ierr);
1014c4762a1bSJed Brown   *min = *max = *mean = 0;
1015c4762a1bSJed Brown   ierr = DMDAGetInfo(da3,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
1016c4762a1bSJed Brown   ierr = DMDAGetCorners(da3,&zs,&ys,&xs,&zm,&ym,&xm);CHKERRQ(ierr);
1017c4762a1bSJed Brown   if (zs != 0 || zm != mz) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Unexpected decomposition");
1018c4762a1bSJed Brown   ierr = DMDAVecGetArray(da3,X3,&x);CHKERRQ(ierr);
1019c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
1020c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
1021c4762a1bSJed Brown       PetscReal u = PetscRealPart(x[i][j][zm-1].u);
1022c4762a1bSJed Brown       RangeUpdate(&umin,&umax,u);
1023c4762a1bSJed Brown       usum += u;
1024c4762a1bSJed Brown     }
1025c4762a1bSJed Brown   }
1026c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da3,X3,&x);CHKERRQ(ierr);
1027c4762a1bSJed Brown   ierr = DMCompositeRestoreAccess(pack,X,&X3,&X2);CHKERRQ(ierr);
1028c4762a1bSJed Brown 
102955b25c41SPierre Jolivet   ierr  = MPI_Allreduce(&umin,min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da3));CHKERRMPI(ierr);
103055b25c41SPierre Jolivet   ierr  = MPI_Allreduce(&umax,max,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da3));CHKERRMPI(ierr);
103155b25c41SPierre Jolivet   ierr  = MPI_Allreduce(&usum,&gusum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da3));CHKERRMPI(ierr);
1032c4762a1bSJed Brown   *mean = PetscRealPart(gusum) / (mx*my);
1033c4762a1bSJed Brown   PetscFunctionReturn(0);
1034c4762a1bSJed Brown }
1035c4762a1bSJed Brown 
1036c4762a1bSJed Brown static PetscErrorCode THISolveStatistics(THI thi,TS ts,PetscInt coarsened,const char name[])
1037c4762a1bSJed Brown {
1038c4762a1bSJed Brown   MPI_Comm       comm;
1039c4762a1bSJed Brown   DM             pack;
1040c4762a1bSJed Brown   Vec            X,X3,X2;
1041c4762a1bSJed Brown   PetscErrorCode ierr;
1042c4762a1bSJed Brown 
1043c4762a1bSJed Brown   PetscFunctionBeginUser;
1044c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)thi,&comm);CHKERRQ(ierr);
1045c4762a1bSJed Brown   ierr = TSGetDM(ts,&pack);CHKERRQ(ierr);
1046c4762a1bSJed Brown   ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
1047c4762a1bSJed Brown   ierr = DMCompositeGetAccess(pack,X,&X3,&X2);CHKERRQ(ierr);
1048c4762a1bSJed Brown   ierr = PetscPrintf(comm,"Solution statistics after solve: %s\n",name);CHKERRQ(ierr);
1049c4762a1bSJed Brown   {
1050c4762a1bSJed Brown     PetscInt            its,lits;
1051c4762a1bSJed Brown     SNESConvergedReason reason;
1052c4762a1bSJed Brown     SNES                snes;
1053c4762a1bSJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1054c4762a1bSJed Brown     ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
1055c4762a1bSJed Brown     ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr);
1056c4762a1bSJed Brown     ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
1057c4762a1bSJed Brown     ierr = PetscPrintf(comm,"%s: Number of SNES iterations = %d, total linear iterations = %d\n",SNESConvergedReasons[reason],its,lits);CHKERRQ(ierr);
1058c4762a1bSJed Brown   }
1059c4762a1bSJed Brown   {
1060c4762a1bSJed Brown     PetscReal   nrm2,tmin[3]={1e100,1e100,1e100},tmax[3]={-1e100,-1e100,-1e100},min[3],max[3];
1061c4762a1bSJed Brown     PetscInt    i,j,m;
1062c4762a1bSJed Brown     PetscScalar *x;
1063c4762a1bSJed Brown     ierr = VecNorm(X3,NORM_2,&nrm2);CHKERRQ(ierr);
1064c4762a1bSJed Brown     ierr = VecGetLocalSize(X3,&m);CHKERRQ(ierr);
1065c4762a1bSJed Brown     ierr = VecGetArray(X3,&x);CHKERRQ(ierr);
1066c4762a1bSJed Brown     for (i=0; i<m; i+=2) {
1067c4762a1bSJed Brown       PetscReal u = PetscRealPart(x[i]),v = PetscRealPart(x[i+1]),c = PetscSqrtReal(u*u+v*v);
1068c4762a1bSJed Brown       tmin[0] = PetscMin(u,tmin[0]);
1069c4762a1bSJed Brown       tmin[1] = PetscMin(v,tmin[1]);
1070c4762a1bSJed Brown       tmin[2] = PetscMin(c,tmin[2]);
1071c4762a1bSJed Brown       tmax[0] = PetscMax(u,tmax[0]);
1072c4762a1bSJed Brown       tmax[1] = PetscMax(v,tmax[1]);
1073c4762a1bSJed Brown       tmax[2] = PetscMax(c,tmax[2]);
1074c4762a1bSJed Brown     }
1075c4762a1bSJed Brown     ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
1076ffc4695bSBarry Smith     ierr = MPI_Allreduce(tmin,min,3,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)thi));CHKERRMPI(ierr);
1077ffc4695bSBarry Smith     ierr = MPI_Allreduce(tmax,max,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)thi));CHKERRMPI(ierr);
1078c4762a1bSJed Brown     /* Dimensionalize to meters/year */
1079c4762a1bSJed Brown     nrm2 *= thi->units->year / thi->units->meter;
1080c4762a1bSJed Brown     for (j=0; j<3; j++) {
1081c4762a1bSJed Brown       min[j] *= thi->units->year / thi->units->meter;
1082c4762a1bSJed Brown       max[j] *= thi->units->year / thi->units->meter;
1083c4762a1bSJed Brown     }
1084c4762a1bSJed Brown     ierr = PetscPrintf(comm,"|X|_2 %g   u in [%g, %g]   v in [%g, %g]   c in [%g, %g] \n",nrm2,min[0],max[0],min[1],max[1],min[2],max[2]);CHKERRQ(ierr);
1085c4762a1bSJed Brown     {
1086c4762a1bSJed Brown       PetscReal umin,umax,umean;
1087c4762a1bSJed Brown       ierr   = THISurfaceStatistics(pack,X,&umin,&umax,&umean);CHKERRQ(ierr);
1088c4762a1bSJed Brown       umin  *= thi->units->year / thi->units->meter;
1089c4762a1bSJed Brown       umax  *= thi->units->year / thi->units->meter;
1090c4762a1bSJed Brown       umean *= thi->units->year / thi->units->meter;
1091c4762a1bSJed Brown       ierr   = PetscPrintf(comm,"Surface statistics: u in [%12.6e, %12.6e] mean %12.6e\n",umin,umax,umean);CHKERRQ(ierr);
1092c4762a1bSJed Brown     }
1093c4762a1bSJed Brown     /* These values stay nondimensional */
1094c4762a1bSJed Brown     ierr = PetscPrintf(comm,"Global eta range   [%g, %g], converged range [%g, %g]\n",thi->eta.min,thi->eta.max,thi->eta.cmin,thi->eta.cmax);CHKERRQ(ierr);
1095c4762a1bSJed Brown     ierr = PetscPrintf(comm,"Global beta2 range [%g, %g], converged range [%g, %g]\n",thi->beta2.min,thi->beta2.max,thi->beta2.cmin,thi->beta2.cmax);CHKERRQ(ierr);
1096c4762a1bSJed Brown   }
1097c4762a1bSJed Brown   ierr = PetscPrintf(comm,"\n");CHKERRQ(ierr);
1098c4762a1bSJed Brown   ierr = DMCompositeRestoreAccess(pack,X,&X3,&X2);CHKERRQ(ierr);
1099c4762a1bSJed Brown   PetscFunctionReturn(0);
1100c4762a1bSJed Brown }
1101c4762a1bSJed Brown 
1102c4762a1bSJed Brown static inline PetscInt DMDALocalIndex3D(DMDALocalInfo *info,PetscInt i,PetscInt j,PetscInt k)
1103c4762a1bSJed Brown {return ((i-info->gzs)*info->gym + (j-info->gys))*info->gxm + (k-info->gxs);}
1104c4762a1bSJed Brown static inline PetscInt DMDALocalIndex2D(DMDALocalInfo *info,PetscInt i,PetscInt j)
1105c4762a1bSJed Brown {return (i-info->gzs)*info->gym + (j-info->gys);}
1106c4762a1bSJed Brown 
1107c4762a1bSJed Brown static PetscErrorCode THIJacobianLocal_Momentum(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,Mat B,Mat Bcpl,THI thi)
1108c4762a1bSJed Brown {
1109c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,zm,i,j,k,q,l,ll;
1110c4762a1bSJed Brown   PetscReal      hx,hy;
1111c4762a1bSJed Brown   PetscErrorCode ierr;
1112c4762a1bSJed Brown 
1113c4762a1bSJed Brown   PetscFunctionBeginUser;
1114c4762a1bSJed Brown   xs = info->zs;
1115c4762a1bSJed Brown   ys = info->ys;
1116c4762a1bSJed Brown   xm = info->zm;
1117c4762a1bSJed Brown   ym = info->ym;
1118c4762a1bSJed Brown   zm = info->xm;
1119c4762a1bSJed Brown   hx = thi->Lx / info->mz;
1120c4762a1bSJed Brown   hy = thi->Ly / info->my;
1121c4762a1bSJed Brown 
1122c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
1123c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
1124c4762a1bSJed Brown       PrmNode pn[4],dpn[4][2];
1125c4762a1bSJed Brown       QuadExtract(prm,i,j,pn);
1126c4762a1bSJed Brown       ierr = QuadComputeGrad4(QuadQDeriv,hx,hy,pn,dpn);CHKERRQ(ierr);
1127c4762a1bSJed Brown       for (k=0; k<zm-1; k++) {
1128c4762a1bSJed Brown         Node        n[8];
1129c4762a1bSJed Brown         PetscReal   zn[8],etabase = 0;
1130c4762a1bSJed Brown         PetscScalar Ke[8*NODE_SIZE][8*NODE_SIZE],Kcpl[8*NODE_SIZE][4*PRMNODE_SIZE];
1131c4762a1bSJed Brown         PetscInt    ls = 0;
1132c4762a1bSJed Brown 
1133c4762a1bSJed Brown         PrmHexGetZ(pn,k,zm,zn);
1134c4762a1bSJed Brown         HexExtract(x,i,j,k,n);
1135c4762a1bSJed Brown         ierr = PetscMemzero(Ke,sizeof(Ke));CHKERRQ(ierr);
1136c4762a1bSJed Brown         ierr = PetscMemzero(Kcpl,sizeof(Kcpl));CHKERRQ(ierr);
1137c4762a1bSJed Brown         if (thi->no_slip && k == 0) {
1138c4762a1bSJed Brown           for (l=0; l<4; l++) n[l].u = n[l].v = 0;
1139c4762a1bSJed Brown           ls = 4;
1140c4762a1bSJed Brown         }
1141c4762a1bSJed Brown         for (q=0; q<8; q++) {
1142c4762a1bSJed Brown           PetscReal   dz[3],phi[8],dphi[8][3],jw,eta,deta;
1143c4762a1bSJed Brown           PetscScalar du[3],dv[3],u,v;
1144c4762a1bSJed Brown           HexGrad(HexQDeriv[q],zn,dz);
1145c4762a1bSJed Brown           HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
1146c4762a1bSJed Brown           PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
1147c4762a1bSJed Brown           jw /= thi->rhog;      /* residuals are scaled by this factor */
1148c4762a1bSJed Brown           if (q == 0) etabase = eta;
1149c4762a1bSJed Brown           for (l=ls; l<8; l++) { /* test functions */
1150c4762a1bSJed Brown             const PetscReal pp=phi[l],*restrict dp = dphi[l];
1151c4762a1bSJed Brown             for (ll=ls; ll<8; ll++) { /* trial functions */
1152c4762a1bSJed Brown               const PetscReal *restrict dpl = dphi[ll];
1153c4762a1bSJed Brown               PetscScalar dgdu,dgdv;
1154c4762a1bSJed 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];
1155c4762a1bSJed 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];
1156c4762a1bSJed Brown               /* Picard part */
1157c4762a1bSJed 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];
1158c4762a1bSJed Brown               Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0];
1159c4762a1bSJed Brown               Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1];
1160c4762a1bSJed 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];
1161c4762a1bSJed Brown               /* extra Newton terms */
1162c4762a1bSJed 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];
1163c4762a1bSJed 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];
1164c4762a1bSJed 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];
1165c4762a1bSJed 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];
1166c4762a1bSJed Brown               /* inertial part */
1167c4762a1bSJed Brown               Ke[l*2+0][ll*2+0] += pp*jw*thi->inertia*pp;
1168c4762a1bSJed Brown               Ke[l*2+1][ll*2+1] += pp*jw*thi->inertia*pp;
1169c4762a1bSJed Brown             }
1170c4762a1bSJed Brown             for (ll=0; ll<4; ll++) { /* Trial functions for surface/bed */
1171c4762a1bSJed Brown               const PetscReal dpl[] = {QuadQDeriv[q%4][ll][0]/hx, QuadQDeriv[q%4][ll][1]/hy}; /* surface = h + b */
1172c4762a1bSJed Brown               Kcpl[FieldIndex(Node,l,u)][FieldIndex(PrmNode,ll,h)] += pp*jw*thi->rhog*dpl[0];
1173c4762a1bSJed Brown               Kcpl[FieldIndex(Node,l,u)][FieldIndex(PrmNode,ll,b)] += pp*jw*thi->rhog*dpl[0];
1174c4762a1bSJed Brown               Kcpl[FieldIndex(Node,l,v)][FieldIndex(PrmNode,ll,h)] += pp*jw*thi->rhog*dpl[1];
1175c4762a1bSJed Brown               Kcpl[FieldIndex(Node,l,v)][FieldIndex(PrmNode,ll,b)] += pp*jw*thi->rhog*dpl[1];
1176c4762a1bSJed Brown             }
1177c4762a1bSJed Brown           }
1178c4762a1bSJed Brown         }
1179c4762a1bSJed Brown         if (k == 0) { /* on a bottom face */
1180c4762a1bSJed Brown           if (thi->no_slip) {
1181c4762a1bSJed Brown             const PetscReal   hz    = PetscRealPart(pn[0].h)/(zm-1);
1182c4762a1bSJed 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);
1183c4762a1bSJed Brown             Ke[0][0] = thi->dirichlet_scale*diagu;
1184c4762a1bSJed Brown             Ke[0][1] = 0;
1185c4762a1bSJed Brown             Ke[1][0] = 0;
1186c4762a1bSJed Brown             Ke[1][1] = thi->dirichlet_scale*diagv;
1187c4762a1bSJed Brown           } else {
1188c4762a1bSJed 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) */
1189c4762a1bSJed Brown               const PetscReal jw = 0.25*hx*hy,*phi = QuadQInterp[q];
1190c4762a1bSJed Brown               PetscScalar     u  =0,v=0,rbeta2=0;
1191c4762a1bSJed Brown               PetscReal       beta2,dbeta2;
1192c4762a1bSJed Brown               for (l=0; l<4; l++) {
1193c4762a1bSJed Brown                 u      += phi[l]*n[l].u;
1194c4762a1bSJed Brown                 v      += phi[l]*n[l].v;
1195c4762a1bSJed Brown                 rbeta2 += phi[l]*pn[l].beta2;
1196c4762a1bSJed Brown               }
1197c4762a1bSJed Brown               THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
1198c4762a1bSJed Brown               for (l=0; l<4; l++) {
1199c4762a1bSJed Brown                 const PetscReal pp = phi[l];
1200c4762a1bSJed Brown                 for (ll=0; ll<4; ll++) {
1201c4762a1bSJed Brown                   const PetscReal ppl = phi[ll];
1202c4762a1bSJed Brown                   Ke[l*2+0][ll*2+0] += pp*jw*beta2*ppl + pp*jw*dbeta2*u*u*ppl;
1203c4762a1bSJed Brown                   Ke[l*2+0][ll*2+1] +=                   pp*jw*dbeta2*u*v*ppl;
1204c4762a1bSJed Brown                   Ke[l*2+1][ll*2+0] +=                   pp*jw*dbeta2*v*u*ppl;
1205c4762a1bSJed Brown                   Ke[l*2+1][ll*2+1] += pp*jw*beta2*ppl + pp*jw*dbeta2*v*v*ppl;
1206c4762a1bSJed Brown                 }
1207c4762a1bSJed Brown               }
1208c4762a1bSJed Brown             }
1209c4762a1bSJed Brown           }
1210c4762a1bSJed Brown         }
1211c4762a1bSJed Brown         {
1212c4762a1bSJed Brown           const PetscInt rc3blocked[8] = {
1213c4762a1bSJed Brown             DMDALocalIndex3D(info,i+0,j+0,k+0),
1214c4762a1bSJed Brown             DMDALocalIndex3D(info,i+1,j+0,k+0),
1215c4762a1bSJed Brown             DMDALocalIndex3D(info,i+1,j+1,k+0),
1216c4762a1bSJed Brown             DMDALocalIndex3D(info,i+0,j+1,k+0),
1217c4762a1bSJed Brown             DMDALocalIndex3D(info,i+0,j+0,k+1),
1218c4762a1bSJed Brown             DMDALocalIndex3D(info,i+1,j+0,k+1),
1219c4762a1bSJed Brown             DMDALocalIndex3D(info,i+1,j+1,k+1),
1220c4762a1bSJed Brown             DMDALocalIndex3D(info,i+0,j+1,k+1)
1221c4762a1bSJed Brown           },col2blocked[PRMNODE_SIZE*4] = {
1222c4762a1bSJed Brown             DMDALocalIndex2D(info,i+0,j+0),
1223c4762a1bSJed Brown             DMDALocalIndex2D(info,i+1,j+0),
1224c4762a1bSJed Brown             DMDALocalIndex2D(info,i+1,j+1),
1225c4762a1bSJed Brown             DMDALocalIndex2D(info,i+0,j+1)
1226c4762a1bSJed Brown           };
1227c4762a1bSJed Brown #if !defined COMPUTE_LOWER_TRIANGULAR /* fill in lower-triangular part, this is really cheap compared to computing the entries */
1228c4762a1bSJed Brown           for (l=0; l<8; l++) {
1229c4762a1bSJed Brown             for (ll=l+1; ll<8; ll++) {
1230c4762a1bSJed Brown               Ke[ll*2+0][l*2+0] = Ke[l*2+0][ll*2+0];
1231c4762a1bSJed Brown               Ke[ll*2+1][l*2+0] = Ke[l*2+0][ll*2+1];
1232c4762a1bSJed Brown               Ke[ll*2+0][l*2+1] = Ke[l*2+1][ll*2+0];
1233c4762a1bSJed Brown               Ke[ll*2+1][l*2+1] = Ke[l*2+1][ll*2+1];
1234c4762a1bSJed Brown             }
1235c4762a1bSJed Brown           }
1236c4762a1bSJed Brown #endif
1237c4762a1bSJed Brown           ierr = MatSetValuesBlockedLocal(B,8,rc3blocked,8,rc3blocked,&Ke[0][0],ADD_VALUES);CHKERRQ(ierr); /* velocity-velocity coupling can use blocked insertion */
1238c4762a1bSJed Brown           {                     /* The off-diagonal part cannot (yet) */
1239c4762a1bSJed Brown             PetscInt row3scalar[NODE_SIZE*8],col2scalar[PRMNODE_SIZE*4];
1240c4762a1bSJed Brown             for (l=0; l<8; l++) for (ll=0; ll<NODE_SIZE; ll++) row3scalar[l*NODE_SIZE+ll] = rc3blocked[l]*NODE_SIZE+ll;
1241c4762a1bSJed Brown             for (l=0; l<4; l++) for (ll=0; ll<PRMNODE_SIZE; ll++) col2scalar[l*PRMNODE_SIZE+ll] = col2blocked[l]*PRMNODE_SIZE+ll;
1242c4762a1bSJed Brown             ierr = MatSetValuesLocal(Bcpl,8*NODE_SIZE,row3scalar,4*PRMNODE_SIZE,col2scalar,&Kcpl[0][0],ADD_VALUES);CHKERRQ(ierr);
1243c4762a1bSJed Brown           }
1244c4762a1bSJed Brown         }
1245c4762a1bSJed Brown       }
1246c4762a1bSJed Brown     }
1247c4762a1bSJed Brown   }
1248c4762a1bSJed Brown   PetscFunctionReturn(0);
1249c4762a1bSJed Brown }
1250c4762a1bSJed Brown 
1251c4762a1bSJed Brown static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info,const Node ***x3,const PrmNode **x2,const PrmNode **xdot2,PetscReal a,Mat B22,Mat B21,THI thi)
1252c4762a1bSJed Brown {
1253c4762a1bSJed Brown   PetscErrorCode ierr;
1254c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,zm,i,j,k;
1255c4762a1bSJed Brown 
1256c4762a1bSJed Brown   PetscFunctionBeginUser;
1257c4762a1bSJed Brown   xs = info->zs;
1258c4762a1bSJed Brown   ys = info->ys;
1259c4762a1bSJed Brown   xm = info->zm;
1260c4762a1bSJed Brown   ym = info->ym;
1261c4762a1bSJed Brown   zm = info->xm;
1262c4762a1bSJed Brown 
1263c4762a1bSJed Brown   if (zm > 1024) SETERRQ(((PetscObject)info->da)->comm,PETSC_ERR_SUP,"Need to allocate more space");
1264c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
1265c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
1266c4762a1bSJed Brown       {                         /* Self-coupling */
1267c4762a1bSJed Brown         const PetscInt    row[]  = {DMDALocalIndex2D(info,i,j)};
1268c4762a1bSJed Brown         const PetscInt    col[]  = {DMDALocalIndex2D(info,i,j)};
1269c4762a1bSJed Brown         const PetscScalar vals[] = {
1270c4762a1bSJed Brown           a,0,0,
1271c4762a1bSJed Brown           0,a,0,
1272c4762a1bSJed Brown           0,0,a
1273c4762a1bSJed Brown         };
1274c4762a1bSJed Brown         ierr = MatSetValuesBlockedLocal(B22,1,row,1,col,vals,INSERT_VALUES);CHKERRQ(ierr);
1275c4762a1bSJed Brown       }
1276c4762a1bSJed Brown       for (k=0; k<zm; k++) {    /* Coupling to velocity problem */
1277c4762a1bSJed Brown         /* Use a cheaper quadrature than for residual evaluation, because it is much sparser */
1278c4762a1bSJed Brown         const PetscInt row[]  = {FieldIndex(PrmNode,DMDALocalIndex2D(info,i,j),h)};
1279c4762a1bSJed Brown         const PetscInt cols[] = {
1280c4762a1bSJed Brown           FieldIndex(Node,DMDALocalIndex3D(info,i-1,j,k),u),
1281c4762a1bSJed Brown           FieldIndex(Node,DMDALocalIndex3D(info,i  ,j,k),u),
1282c4762a1bSJed Brown           FieldIndex(Node,DMDALocalIndex3D(info,i+1,j,k),u),
1283c4762a1bSJed Brown           FieldIndex(Node,DMDALocalIndex3D(info,i,j-1,k),v),
1284c4762a1bSJed Brown           FieldIndex(Node,DMDALocalIndex3D(info,i,j  ,k),v),
1285c4762a1bSJed Brown           FieldIndex(Node,DMDALocalIndex3D(info,i,j+1,k),v)
1286c4762a1bSJed Brown         };
1287c4762a1bSJed Brown         const PetscScalar
1288c4762a1bSJed Brown           w  = (k && k<zm-1) ? 0.5 : 0.25,
1289c4762a1bSJed Brown           hW = w*(x2[i-1][j  ].h+x2[i  ][j  ].h)/(zm-1.),
1290c4762a1bSJed Brown           hE = w*(x2[i  ][j  ].h+x2[i+1][j  ].h)/(zm-1.),
1291c4762a1bSJed Brown           hS = w*(x2[i  ][j-1].h+x2[i  ][j  ].h)/(zm-1.),
1292c4762a1bSJed Brown           hN = w*(x2[i  ][j  ].h+x2[i  ][j+1].h)/(zm-1.);
1293c4762a1bSJed Brown         PetscScalar *vals,
1294c4762a1bSJed Brown                      vals_upwind[] = {((PetscRealPart(x3[i][j][k].u) > 0) ? -hW : 0),
1295c4762a1bSJed Brown                                       ((PetscRealPart(x3[i][j][k].u) > 0) ? +hE : -hW),
1296c4762a1bSJed Brown                                       ((PetscRealPart(x3[i][j][k].u) > 0) ?  0  : +hE),
1297c4762a1bSJed Brown                                       ((PetscRealPart(x3[i][j][k].v) > 0) ? -hS : 0),
1298c4762a1bSJed Brown                                       ((PetscRealPart(x3[i][j][k].v) > 0) ? +hN : -hS),
1299c4762a1bSJed Brown                                       ((PetscRealPart(x3[i][j][k].v) > 0) ?  0  : +hN)},
1300c4762a1bSJed Brown                      vals_centered[] = {-0.5*hW, 0.5*(-hW+hE), 0.5*hE,
1301c4762a1bSJed Brown                                         -0.5*hS, 0.5*(-hS+hN), 0.5*hN};
1302c4762a1bSJed Brown         vals = 1 ? vals_upwind : vals_centered;
1303c4762a1bSJed Brown         if (k == 0) {
1304c4762a1bSJed Brown           Node derate;
1305c4762a1bSJed Brown           THIErosion(thi,&x3[i][j][0],NULL,&derate);
1306c4762a1bSJed Brown           vals[1] -= derate.u;
1307c4762a1bSJed Brown           vals[4] -= derate.v;
1308c4762a1bSJed Brown         }
1309c4762a1bSJed Brown         ierr = MatSetValuesLocal(B21,1,row,6,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1310c4762a1bSJed Brown       }
1311c4762a1bSJed Brown     }
1312c4762a1bSJed Brown   }
1313c4762a1bSJed Brown   PetscFunctionReturn(0);
1314c4762a1bSJed Brown }
1315c4762a1bSJed Brown 
1316c4762a1bSJed Brown static PetscErrorCode THIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
1317c4762a1bSJed Brown {
1318c4762a1bSJed Brown   PetscErrorCode ierr;
1319c4762a1bSJed Brown   THI            thi = (THI)ctx;
1320c4762a1bSJed Brown   DM             pack,da3,da2;
1321c4762a1bSJed Brown   Vec            X3,X2,Xdot2;
1322c4762a1bSJed Brown   Mat            B11,B12,B21,B22;
1323c4762a1bSJed Brown   DMDALocalInfo  info3;
1324c4762a1bSJed Brown   IS             *isloc;
1325c4762a1bSJed Brown   const Node     ***x3;
1326c4762a1bSJed Brown   const PrmNode  **x2,**xdot2;
1327c4762a1bSJed Brown 
1328c4762a1bSJed Brown   PetscFunctionBeginUser;
1329c4762a1bSJed Brown   ierr = TSGetDM(ts,&pack);CHKERRQ(ierr);
1330c4762a1bSJed Brown   ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr);
1331c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da3,&info3);CHKERRQ(ierr);
1332c4762a1bSJed Brown   ierr = DMCompositeGetLocalVectors(pack,&X3,&X2);CHKERRQ(ierr);
1333c4762a1bSJed Brown   ierr = DMCompositeGetLocalVectors(pack,NULL,&Xdot2);CHKERRQ(ierr);
1334c4762a1bSJed Brown   ierr = DMCompositeScatter(pack,X,X3,X2);CHKERRQ(ierr);
1335c4762a1bSJed Brown   ierr = THIFixGhosts(thi,da3,da2,X3,X2);CHKERRQ(ierr);
1336c4762a1bSJed Brown   ierr = DMCompositeScatter(pack,Xdot,NULL,Xdot2);CHKERRQ(ierr);
1337c4762a1bSJed Brown 
1338c4762a1bSJed Brown   ierr = MatZeroEntries(B);CHKERRQ(ierr);
1339c4762a1bSJed Brown 
1340c4762a1bSJed Brown   ierr = DMCompositeGetLocalISs(pack,&isloc);CHKERRQ(ierr);
1341c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(B,isloc[0],isloc[0],&B11);CHKERRQ(ierr);
1342c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(B,isloc[0],isloc[1],&B12);CHKERRQ(ierr);
1343c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(B,isloc[1],isloc[0],&B21);CHKERRQ(ierr);
1344c4762a1bSJed Brown   ierr = MatGetLocalSubMatrix(B,isloc[1],isloc[1],&B22);CHKERRQ(ierr);
1345c4762a1bSJed Brown 
1346c4762a1bSJed Brown   ierr = DMDAVecGetArray(da3,X3,&x3);CHKERRQ(ierr);
1347c4762a1bSJed Brown   ierr = DMDAVecGetArray(da2,X2,&x2);CHKERRQ(ierr);
1348c4762a1bSJed Brown   ierr = DMDAVecGetArray(da2,Xdot2,&xdot2);CHKERRQ(ierr);
1349c4762a1bSJed Brown 
1350c4762a1bSJed Brown   ierr = THIJacobianLocal_Momentum(&info3,x3,x2,B11,B12,thi);CHKERRQ(ierr);
1351c4762a1bSJed Brown 
1352c4762a1bSJed Brown   /* Need to switch from ADD_VALUES to INSERT_VALUES */
1353c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
1354c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
1355c4762a1bSJed Brown 
1356c4762a1bSJed Brown   ierr = THIJacobianLocal_2D(&info3,x3,x2,xdot2,a,B22,B21,thi);CHKERRQ(ierr);
1357c4762a1bSJed Brown 
1358c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da3,X3,&x3);CHKERRQ(ierr);
1359c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da2,X2,&x2);CHKERRQ(ierr);
1360c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da2,Xdot2,&xdot2);CHKERRQ(ierr);
1361c4762a1bSJed Brown 
1362c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(B,isloc[0],isloc[0],&B11);CHKERRQ(ierr);
1363c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(B,isloc[0],isloc[1],&B12);CHKERRQ(ierr);
1364c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(B,isloc[1],isloc[0],&B21);CHKERRQ(ierr);
1365c4762a1bSJed Brown   ierr = MatRestoreLocalSubMatrix(B,isloc[1],isloc[1],&B22);CHKERRQ(ierr);
1366c4762a1bSJed Brown   ierr = ISDestroy(&isloc[0]);CHKERRQ(ierr);
1367c4762a1bSJed Brown   ierr = ISDestroy(&isloc[1]);CHKERRQ(ierr);
1368c4762a1bSJed Brown   ierr = PetscFree(isloc);CHKERRQ(ierr);
1369c4762a1bSJed Brown 
1370c4762a1bSJed Brown   ierr = DMCompositeRestoreLocalVectors(pack,&X3,&X2);CHKERRQ(ierr);
1371c4762a1bSJed Brown   ierr = DMCompositeRestoreLocalVectors(pack,0,&Xdot2);CHKERRQ(ierr);
1372c4762a1bSJed Brown 
1373c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1374c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1375c4762a1bSJed Brown   if (A != B) {
1376c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1377c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1378c4762a1bSJed Brown   }
1379c4762a1bSJed Brown   if (thi->verbose) {ierr = THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
1380c4762a1bSJed Brown   PetscFunctionReturn(0);
1381c4762a1bSJed Brown }
1382c4762a1bSJed Brown 
1383c4762a1bSJed Brown /* VTK's XML formats are so brain-dead that they can't handle multiple grids in the same file.  Since the communication
1384c4762a1bSJed Brown  * can be shared between the two grids, we write two files at once, one for velocity (living on a 3D grid defined by
1385c4762a1bSJed Brown  * h=thickness and b=bed) and another for all properties living on the 2D grid.
1386c4762a1bSJed Brown  */
1387c4762a1bSJed Brown static PetscErrorCode THIDAVecView_VTK_XML(THI thi,DM pack,Vec X,const char filename[],const char filename2[])
1388c4762a1bSJed Brown {
1389c4762a1bSJed Brown   const PetscInt dof   = NODE_SIZE,dof2 = PRMNODE_SIZE;
1390c4762a1bSJed Brown   Units          units = thi->units;
1391c4762a1bSJed Brown   MPI_Comm       comm;
1392c4762a1bSJed Brown   PetscErrorCode ierr;
1393c4762a1bSJed Brown   PetscViewer    viewer3,viewer2;
1394c4762a1bSJed Brown   PetscMPIInt    rank,size,tag,nn,nmax,nn2,nmax2;
1395c4762a1bSJed Brown   PetscInt       mx,my,mz,r,range[6];
1396c4762a1bSJed Brown   PetscScalar    *x,*x2;
1397c4762a1bSJed Brown   DM             da3,da2;
1398c4762a1bSJed Brown   Vec            X3,X2;
1399c4762a1bSJed Brown 
1400c4762a1bSJed Brown   PetscFunctionBeginUser;
1401c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)thi,&comm);CHKERRQ(ierr);
1402c4762a1bSJed Brown   ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr);
1403c4762a1bSJed Brown   ierr = DMCompositeGetAccess(pack,X,&X3,&X2);CHKERRQ(ierr);
1404c4762a1bSJed Brown   ierr = DMDAGetInfo(da3,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
1405ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
1406ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
1407c4762a1bSJed Brown   ierr = PetscViewerASCIIOpen(comm,filename,&viewer3);CHKERRQ(ierr);
1408c4762a1bSJed Brown   ierr = PetscViewerASCIIOpen(comm,filename2,&viewer2);CHKERRQ(ierr);
1409c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer3,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
1410c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer2,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
1411c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer3,"  <StructuredGrid WholeExtent=\"%d %d %d %d %d %d\">\n",0,mz-1,0,my-1,0,mx-1);CHKERRQ(ierr);
1412c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer2,"  <StructuredGrid WholeExtent=\"%d %d %d %d %d %d\">\n",0,0,0,my-1,0,mx-1);CHKERRQ(ierr);
1413c4762a1bSJed Brown 
1414c4762a1bSJed Brown   ierr = DMDAGetCorners(da3,range,range+1,range+2,range+3,range+4,range+5);CHKERRQ(ierr);
1415c4762a1bSJed Brown   ierr = PetscMPIIntCast(range[3]*range[4]*range[5]*dof,&nn);CHKERRQ(ierr);
1416ffc4695bSBarry Smith   ierr = MPI_Reduce(&nn,&nmax,1,MPI_INT,MPI_MAX,0,comm);CHKERRMPI(ierr);
1417c4762a1bSJed Brown   ierr = PetscMPIIntCast(range[4]*range[5]*dof2,&nn2);CHKERRQ(ierr);
1418ffc4695bSBarry Smith   ierr = MPI_Reduce(&nn2,&nmax2,1,MPI_INT,MPI_MAX,0,comm);CHKERRMPI(ierr);
1419c4762a1bSJed Brown   tag  = ((PetscObject)viewer3)->tag;
1420c4762a1bSJed Brown   ierr = VecGetArrayRead(X3,(const PetscScalar**)&x);CHKERRQ(ierr);
1421c4762a1bSJed Brown   ierr = VecGetArrayRead(X2,(const PetscScalar**)&x2);CHKERRQ(ierr);
1422dd400576SPatrick Sanan   if (rank == 0) {
1423c4762a1bSJed Brown     PetscScalar *array,*array2;
1424c4762a1bSJed Brown     ierr = PetscMalloc2(nmax,&array,nmax2,&array2);CHKERRQ(ierr);
1425c4762a1bSJed Brown     for (r=0; r<size; r++) {
1426c4762a1bSJed Brown       PetscInt    i,j,k,f,xs,xm,ys,ym,zs,zm;
1427c4762a1bSJed Brown       Node        *y3;
1428c4762a1bSJed Brown       PetscScalar (*y2)[PRMNODE_SIZE];
1429c4762a1bSJed Brown       MPI_Status status;
1430c4762a1bSJed Brown 
1431c4762a1bSJed Brown       if (r) {
1432ffc4695bSBarry Smith         ierr = MPI_Recv(range,6,MPIU_INT,r,tag,comm,MPI_STATUS_IGNORE);CHKERRMPI(ierr);
1433c4762a1bSJed Brown       }
1434c4762a1bSJed Brown       zs = range[0];ys = range[1];xs = range[2];zm = range[3];ym = range[4];xm = range[5];
1435c4762a1bSJed Brown       if (xm*ym*zm*dof > nmax) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not happen");
1436c4762a1bSJed Brown       if (r) {
1437ffc4695bSBarry Smith         ierr = MPI_Recv(array,nmax,MPIU_SCALAR,r,tag,comm,&status);CHKERRMPI(ierr);
1438ffc4695bSBarry Smith         ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRMPI(ierr);
1439c4762a1bSJed Brown         if (nn != xm*ym*zm*dof) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"corrupt da3 send");
1440c4762a1bSJed Brown         y3   = (Node*)array;
1441ffc4695bSBarry Smith         ierr = MPI_Recv(array2,nmax2,MPIU_SCALAR,r,tag,comm,&status);CHKERRMPI(ierr);
1442ffc4695bSBarry Smith         ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn2);CHKERRMPI(ierr);
1443c4762a1bSJed Brown         if (nn2 != xm*ym*dof2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"corrupt da2 send");
1444c4762a1bSJed Brown         y2 = (PetscScalar(*)[PRMNODE_SIZE])array2;
1445c4762a1bSJed Brown       } else {
1446c4762a1bSJed Brown         y3 = (Node*)x;
1447c4762a1bSJed Brown         y2 = (PetscScalar(*)[PRMNODE_SIZE])x2;
1448c4762a1bSJed Brown       }
1449c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer3,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",zs,zs+zm-1,ys,ys+ym-1,xs,xs+xm-1);CHKERRQ(ierr);
1450c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer2,"    <Piece Extent=\"%d %d %D %D %D %D\">\n",0,0,ys,ys+ym-1,xs,xs+xm-1);CHKERRQ(ierr);
1451c4762a1bSJed Brown 
1452c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer3,"      <Points>\n");CHKERRQ(ierr);
1453c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer2,"      <Points>\n");CHKERRQ(ierr);
1454c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer3,"        <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n");CHKERRQ(ierr);
1455c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer2,"        <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n");CHKERRQ(ierr);
1456c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
1457c4762a1bSJed Brown         for (j=ys; j<ys+ym; j++) {
1458c4762a1bSJed Brown           PetscReal
1459c4762a1bSJed Brown             xx = thi->Lx*i/mx,
1460c4762a1bSJed Brown             yy = thi->Ly*j/my,
1461c4762a1bSJed Brown             b  = PetscRealPart(y2[i*ym+j][FieldOffset(PrmNode,b)]),
1462c4762a1bSJed Brown             h  = PetscRealPart(y2[i*ym+j][FieldOffset(PrmNode,h)]);
1463c4762a1bSJed Brown           for (k=zs; k<zs+zm; k++) {
1464c4762a1bSJed Brown             PetscReal zz = b + h*k/(mz-1.);
1465c4762a1bSJed Brown             ierr = PetscViewerASCIIPrintf(viewer3,"%f %f %f\n",xx,yy,zz);CHKERRQ(ierr);
1466c4762a1bSJed Brown           }
1467c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer2,"%f %f %f\n",xx,yy,(double)0.0);CHKERRQ(ierr);
1468c4762a1bSJed Brown         }
1469c4762a1bSJed Brown       }
1470c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer3,"        </DataArray>\n");CHKERRQ(ierr);
1471c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer2,"        </DataArray>\n");CHKERRQ(ierr);
1472c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer3,"      </Points>\n");CHKERRQ(ierr);
1473c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer2,"      </Points>\n");CHKERRQ(ierr);
1474c4762a1bSJed Brown 
1475c4762a1bSJed Brown       {                         /* Velocity and rank (3D) */
1476c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer3,"      <PointData>\n");CHKERRQ(ierr);
1477c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer3,"        <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n");CHKERRQ(ierr);
1478c4762a1bSJed Brown         for (i=0; i<nn/dof; i++) {
1479c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer3,"%f %f %f\n",PetscRealPart(y3[i].u)*units->year/units->meter,PetscRealPart(y3[i].v)*units->year/units->meter,0.0);CHKERRQ(ierr);
1480c4762a1bSJed Brown         }
1481c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer3,"        </DataArray>\n");CHKERRQ(ierr);
1482c4762a1bSJed Brown 
1483c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer3,"        <DataArray type=\"Int32\" Name=\"rank\" NumberOfComponents=\"1\" format=\"ascii\">\n");CHKERRQ(ierr);
1484c4762a1bSJed Brown         for (i=0; i<nn; i+=dof) {
1485c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer3,"%D\n",r);CHKERRQ(ierr);
1486c4762a1bSJed Brown         }
1487c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer3,"        </DataArray>\n");CHKERRQ(ierr);
1488c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer3,"      </PointData>\n");CHKERRQ(ierr);
1489c4762a1bSJed Brown       }
1490c4762a1bSJed Brown 
1491c4762a1bSJed Brown       {                         /* 2D */
1492c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer2,"      <PointData>\n");CHKERRQ(ierr);
1493c4762a1bSJed Brown         for (f=0; f<PRMNODE_SIZE; f++) {
1494c4762a1bSJed Brown           const char *fieldname;
1495c4762a1bSJed Brown           ierr = DMDAGetFieldName(da2,f,&fieldname);CHKERRQ(ierr);
1496c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer2,"        <DataArray type=\"Float32\" Name=\"%s\" format=\"ascii\">\n",fieldname);CHKERRQ(ierr);
1497c4762a1bSJed Brown           for (i=0; i<nn2/PRMNODE_SIZE; i++) {
1498c4762a1bSJed Brown             ierr = PetscViewerASCIIPrintf(viewer2,"%g\n",y2[i][f]);CHKERRQ(ierr);
1499c4762a1bSJed Brown           }
1500c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer2,"        </DataArray>\n");CHKERRQ(ierr);
1501c4762a1bSJed Brown         }
1502c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer2,"      </PointData>\n");CHKERRQ(ierr);
1503c4762a1bSJed Brown       }
1504c4762a1bSJed Brown 
1505c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer3,"    </Piece>\n");CHKERRQ(ierr);
1506c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer2,"    </Piece>\n");CHKERRQ(ierr);
1507c4762a1bSJed Brown     }
1508c4762a1bSJed Brown     ierr = PetscFree2(array,array2);CHKERRQ(ierr);
1509c4762a1bSJed Brown   } else {
1510ffc4695bSBarry Smith     ierr = MPI_Send(range,6,MPIU_INT,0,tag,comm);CHKERRMPI(ierr);
1511ffc4695bSBarry Smith     ierr = MPI_Send(x,nn,MPIU_SCALAR,0,tag,comm);CHKERRMPI(ierr);
1512ffc4695bSBarry Smith     ierr = MPI_Send(x2,nn2,MPIU_SCALAR,0,tag,comm);CHKERRMPI(ierr);
1513c4762a1bSJed Brown   }
1514c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X3,(const PetscScalar**)&x);CHKERRQ(ierr);
1515c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X2,(const PetscScalar**)&x2);CHKERRQ(ierr);
1516c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer3,"  </StructuredGrid>\n");CHKERRQ(ierr);
1517c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer2,"  </StructuredGrid>\n");CHKERRQ(ierr);
1518c4762a1bSJed Brown 
1519c4762a1bSJed Brown   ierr = DMCompositeRestoreAccess(pack,X,&X3,&X2);CHKERRQ(ierr);
1520c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer3,"</VTKFile>\n");CHKERRQ(ierr);
1521c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(viewer2,"</VTKFile>\n");CHKERRQ(ierr);
1522c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer3);CHKERRQ(ierr);
1523c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer2);CHKERRQ(ierr);
1524c4762a1bSJed Brown   PetscFunctionReturn(0);
1525c4762a1bSJed Brown }
1526c4762a1bSJed Brown 
1527c4762a1bSJed Brown static PetscErrorCode THITSMonitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
1528c4762a1bSJed Brown {
1529c4762a1bSJed Brown   PetscErrorCode ierr;
1530c4762a1bSJed Brown   THI            thi = (THI)ctx;
1531c4762a1bSJed Brown   DM             pack;
1532c4762a1bSJed Brown   char           filename3[PETSC_MAX_PATH_LEN],filename2[PETSC_MAX_PATH_LEN];
1533c4762a1bSJed Brown 
1534c4762a1bSJed Brown   PetscFunctionBeginUser;
1535c4762a1bSJed Brown   if (step < 0) PetscFunctionReturn(0); /* negative one is used to indicate an interpolated solution */
1536c4762a1bSJed Brown   ierr = PetscPrintf(PetscObjectComm((PetscObject)ts),"%3D: t=%g\n",step,(double)t);CHKERRQ(ierr);
1537c4762a1bSJed Brown   if (thi->monitor_interval && step % thi->monitor_interval) PetscFunctionReturn(0);
1538c4762a1bSJed Brown   ierr = TSGetDM(ts,&pack);CHKERRQ(ierr);
1539c4762a1bSJed Brown   ierr = PetscSNPrintf(filename3,sizeof(filename3),"%s-3d-%03d.vts",thi->monitor_basename,step);CHKERRQ(ierr);
1540c4762a1bSJed Brown   ierr = PetscSNPrintf(filename2,sizeof(filename2),"%s-2d-%03d.vts",thi->monitor_basename,step);CHKERRQ(ierr);
1541c4762a1bSJed Brown   ierr = THIDAVecView_VTK_XML(thi,pack,X,filename3,filename2);CHKERRQ(ierr);
1542c4762a1bSJed Brown   PetscFunctionReturn(0);
1543c4762a1bSJed Brown }
1544c4762a1bSJed Brown 
1545c4762a1bSJed Brown static PetscErrorCode THICreateDM3d(THI thi,DM *dm3d)
1546c4762a1bSJed Brown {
1547c4762a1bSJed Brown   MPI_Comm       comm;
1548c4762a1bSJed Brown   PetscInt       M    = 3,N = 3,P = 2;
1549c4762a1bSJed Brown   DM             da;
1550c4762a1bSJed Brown   PetscErrorCode ierr;
1551c4762a1bSJed Brown 
1552c4762a1bSJed Brown   PetscFunctionBeginUser;
1553c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)thi,&comm);CHKERRQ(ierr);
1554c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"Grid resolution options","");CHKERRQ(ierr);
1555c4762a1bSJed Brown   {
1556c4762a1bSJed Brown     ierr = PetscOptionsInt("-M","Number of elements in x-direction on coarse level","",M,&M,NULL);CHKERRQ(ierr);
1557c4762a1bSJed Brown     N    = M;
1558c4762a1bSJed Brown     ierr = PetscOptionsInt("-N","Number of elements in y-direction on coarse level (if different from M)","",N,&N,NULL);CHKERRQ(ierr);
1559c4762a1bSJed Brown     ierr = PetscOptionsInt("-P","Number of elements in z-direction on coarse level","",P,&P,NULL);CHKERRQ(ierr);
1560c4762a1bSJed Brown   }
1561c4762a1bSJed Brown   ierr  = PetscOptionsEnd();CHKERRQ(ierr);
1562c4762a1bSJed Brown   ierr  = 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);CHKERRQ(ierr);
1563c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
1564c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
1565c4762a1bSJed Brown   ierr  = DMDASetFieldName(da,0,"x-velocity");CHKERRQ(ierr);
1566c4762a1bSJed Brown   ierr  = DMDASetFieldName(da,1,"y-velocity");CHKERRQ(ierr);
1567c4762a1bSJed Brown   *dm3d = da;
1568c4762a1bSJed Brown   PetscFunctionReturn(0);
1569c4762a1bSJed Brown }
1570c4762a1bSJed Brown 
1571c4762a1bSJed Brown int main(int argc,char *argv[])
1572c4762a1bSJed Brown {
1573c4762a1bSJed Brown   MPI_Comm       comm;
1574c4762a1bSJed Brown   DM             pack,da3,da2;
1575c4762a1bSJed Brown   TS             ts;
1576c4762a1bSJed Brown   THI            thi;
1577c4762a1bSJed Brown   Vec            X;
1578c4762a1bSJed Brown   Mat            B;
1579c4762a1bSJed Brown   PetscInt       i,steps;
1580c4762a1bSJed Brown   PetscReal      ftime;
1581c4762a1bSJed Brown   PetscErrorCode ierr;
1582c4762a1bSJed Brown 
1583c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
1584c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1585c4762a1bSJed Brown 
1586c4762a1bSJed Brown   ierr = THICreate(comm,&thi);CHKERRQ(ierr);
1587c4762a1bSJed Brown   ierr = THICreateDM3d(thi,&da3);CHKERRQ(ierr);
1588c4762a1bSJed Brown   {
1589c4762a1bSJed Brown     PetscInt        Mx,My,mx,my,s;
1590c4762a1bSJed Brown     DMDAStencilType st;
1591c4762a1bSJed Brown     ierr = DMDAGetInfo(da3,0, 0,&My,&Mx, 0,&my,&mx, 0,&s,0,0,0,&st);CHKERRQ(ierr);
1592c4762a1bSJed Brown     ierr = DMDACreate2d(PetscObjectComm((PetscObject)thi),DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,My,Mx,my,mx,sizeof(PrmNode)/sizeof(PetscScalar),s,0,0,&da2);CHKERRQ(ierr);
1593c4762a1bSJed Brown     ierr = DMSetUp(da2);CHKERRQ(ierr);
1594c4762a1bSJed Brown   }
1595c4762a1bSJed Brown 
1596c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)da3,"3D_Velocity");CHKERRQ(ierr);
1597c4762a1bSJed Brown   ierr = DMSetOptionsPrefix(da3,"f3d_");CHKERRQ(ierr);
1598c4762a1bSJed Brown   ierr = DMDASetFieldName(da3,0,"u");CHKERRQ(ierr);
1599c4762a1bSJed Brown   ierr = DMDASetFieldName(da3,1,"v");CHKERRQ(ierr);
1600c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)da2,"2D_Fields");CHKERRQ(ierr);
1601c4762a1bSJed Brown   ierr = DMSetOptionsPrefix(da2,"f2d_");CHKERRQ(ierr);
1602c4762a1bSJed Brown   ierr = DMDASetFieldName(da2,0,"b");CHKERRQ(ierr);
1603c4762a1bSJed Brown   ierr = DMDASetFieldName(da2,1,"h");CHKERRQ(ierr);
1604c4762a1bSJed Brown   ierr = DMDASetFieldName(da2,2,"beta2");CHKERRQ(ierr);
1605c4762a1bSJed Brown   ierr = DMCompositeCreate(comm,&pack);CHKERRQ(ierr);
1606c4762a1bSJed Brown   ierr = DMCompositeAddDM(pack,da3);CHKERRQ(ierr);
1607c4762a1bSJed Brown   ierr = DMCompositeAddDM(pack,da2);CHKERRQ(ierr);
1608c4762a1bSJed Brown   ierr = DMDestroy(&da3);CHKERRQ(ierr);
1609c4762a1bSJed Brown   ierr = DMDestroy(&da2);CHKERRQ(ierr);
1610c4762a1bSJed Brown   ierr = DMSetUp(pack);CHKERRQ(ierr);
1611c4762a1bSJed Brown   ierr = DMCreateMatrix(pack,&B);CHKERRQ(ierr);
1612c4762a1bSJed Brown   ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
1613c4762a1bSJed Brown   ierr = MatSetOptionsPrefix(B,"thi_");CHKERRQ(ierr);
1614c4762a1bSJed Brown 
1615c4762a1bSJed Brown   for (i=0; i<thi->nlevels; i++) {
1616c4762a1bSJed Brown     PetscReal Lx = thi->Lx / thi->units->meter,Ly = thi->Ly / thi->units->meter,Lz = thi->Lz / thi->units->meter;
1617c4762a1bSJed Brown     PetscInt  Mx,My,Mz;
1618c4762a1bSJed Brown     ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr);
1619c4762a1bSJed Brown     ierr = DMDAGetInfo(da3,0, &Mz,&My,&Mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
1620c4762a1bSJed Brown     ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Level %D 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));CHKERRQ(ierr);
1621c4762a1bSJed Brown   }
1622c4762a1bSJed Brown 
1623c4762a1bSJed Brown   ierr = DMCreateGlobalVector(pack,&X);CHKERRQ(ierr);
1624c4762a1bSJed Brown   ierr = THIInitial(thi,pack,X);CHKERRQ(ierr);
1625c4762a1bSJed Brown 
1626c4762a1bSJed Brown   ierr = TSCreate(comm,&ts);CHKERRQ(ierr);
1627c4762a1bSJed Brown   ierr = TSSetDM(ts,pack);CHKERRQ(ierr);
1628c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
1629c4762a1bSJed Brown   ierr = TSMonitorSet(ts,THITSMonitor,thi,NULL);CHKERRQ(ierr);
1630c4762a1bSJed Brown   ierr = TSSetType(ts,TSTHETA);CHKERRQ(ierr);
1631c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,THIFunction,thi);CHKERRQ(ierr);
1632c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,B,B,THIJacobian,thi);CHKERRQ(ierr);
1633c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,10.0);CHKERRQ(ierr);
1634c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
1635c4762a1bSJed Brown   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
1636c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,1e-3);CHKERRQ(ierr);
1637c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
1638c4762a1bSJed Brown 
1639c4762a1bSJed Brown   ierr = TSSolve(ts,X);CHKERRQ(ierr);
1640c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
1641c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
1642c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Steps %D  final time %g\n",steps,(double)ftime);CHKERRQ(ierr);
1643c4762a1bSJed Brown 
1644c4762a1bSJed Brown   if (0) {ierr = THISolveStatistics(thi,ts,0,"Full");CHKERRQ(ierr);}
1645c4762a1bSJed Brown 
1646c4762a1bSJed Brown   {
1647c4762a1bSJed Brown     PetscBool flg;
1648c4762a1bSJed Brown     char      filename[PETSC_MAX_PATH_LEN] = "";
1649c4762a1bSJed Brown     ierr = PetscOptionsGetString(NULL,NULL,"-o",filename,sizeof(filename),&flg);CHKERRQ(ierr);
1650c4762a1bSJed Brown     if (flg) {
1651c4762a1bSJed Brown       ierr = THIDAVecView_VTK_XML(thi,pack,X,filename,NULL);CHKERRQ(ierr);
1652c4762a1bSJed Brown     }
1653c4762a1bSJed Brown   }
1654c4762a1bSJed Brown 
1655c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
1656c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
1657c4762a1bSJed Brown   ierr = DMDestroy(&pack);CHKERRQ(ierr);
1658c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
1659c4762a1bSJed Brown   ierr = THIDestroy(&thi);CHKERRQ(ierr);
1660c4762a1bSJed Brown   ierr = PetscFinalize();
1661c4762a1bSJed Brown   return ierr;
1662c4762a1bSJed Brown }
1663