Actual source code: ex48.c
petsc-3.7.3 2016-08-01
1: static const char help[] = "Toy hydrostatic ice flow with multigrid in 3D.\n\
2: \n\
3: Solves the hydrostatic (aka Blatter/Pattyn/First Order) equations for ice sheet flow\n\
4: using multigrid. The ice uses a power-law rheology with \"Glen\" exponent 3 (corresponds\n\
5: to p=4/3 in a p-Laplacian). The focus is on ISMIP-HOM experiments which assume periodic\n\
6: boundary conditions in the x- and y-directions.\n\
7: \n\
8: Equations are rescaled so that the domain size and solution are O(1), details of this scaling\n\
9: can be controlled by the options -units_meter, -units_second, and -units_kilogram.\n\
10: \n\
11: A VTK StructuredGrid output file can be written using the option -o filename.vts\n\
12: \n\n";
14: /*
15: The equations for horizontal velocity (u,v) are
17: - [eta (4 u_x + 2 v_y)]_x - [eta (u_y + v_x)]_y - [eta u_z]_z + rho g s_x = 0
18: - [eta (4 v_y + 2 u_x)]_y - [eta (u_y + v_x)]_x - [eta v_z]_z + rho g s_y = 0
20: where
22: eta = B/2 (epsilon + gamma)^((p-2)/2)
24: is the nonlinear effective viscosity with regularization epsilon and hardness parameter B,
25: written in terms of the second invariant
27: 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
29: The surface boundary conditions are the natural conditions. The basal boundary conditions
30: are either no-slip, or Navier (linear) slip with spatially variant friction coefficient beta^2.
32: In the code, the equations for (u,v) are multiplied through by 1/(rho g) so that residuals are O(1).
34: The discretization is Q1 finite elements, managed by a DMDA. The grid is never distorted in the
35: map (x,y) plane, but the bed and surface may be bumpy. This is handled as usual in FEM, through
36: the Jacobian of the coordinate transformation from a reference element to the physical element.
38: Since ice-flow is tightly coupled in the z-direction (within columns), the DMDA is managed
39: specially so that columns are never distributed, and are always contiguous in memory.
40: This amounts to reversing the meaning of X,Y,Z compared to the DMDA's internal interpretation,
41: and then indexing as vec[i][j][k]. The exotic coarse spaces require 2D DMDAs which are made to
42: use compatible domain decomposition relative to the 3D DMDAs.
44: There are two compile-time options:
46: NO_SSE2:
47: If the host supports SSE2, we use integration code that has been vectorized with SSE2
48: intrinsics, unless this macro is defined. The intrinsics speed up integration by about
49: 30% on my architecture (P8700, gcc-4.5 snapshot).
51: COMPUTE_LOWER_TRIANGULAR:
52: The element matrices we assemble are lower-triangular so it is not necessary to compute
53: all entries explicitly. If this macro is defined, the lower-triangular entries are
54: computed explicitly.
56: */
58: #if defined(PETSC_APPLE_FRAMEWORK)
59: #import <PETSc/petscsnes.h>
60: #import <PETSc/petsc/private/dmdaimpl.h> /* There is not yet a public interface to manipulate dm->ops */
61: #else
62: #include <petscsnes.h>
63: #include <petsc/private/dmdaimpl.h> /* There is not yet a public interface to manipulate dm->ops */
64: #endif
65: #include <ctype.h> /* toupper() */
67: #if defined(__cplusplus)
68: /* c++ cannot handle [_restrict_] notation like C does */
69: #undef PETSC_RESTRICT
70: #define PETSC_RESTRICT
71: #endif
73: #if defined __SSE2__
74: # include <emmintrin.h>
75: #endif
77: /* The SSE2 kernels are only for PetscScalar=double on architectures that support it */
78: #define USE_SSE2_KERNELS (!defined NO_SSE2 \
79: && !defined PETSC_USE_COMPLEX \
80: && !defined PETSC_USE_REAL_SINGLE \
81: && !defined PETSC_USE_REAL___FLOAT128 \
82: && defined __SSE2__)
84: static PetscClassId THI_CLASSID;
86: typedef enum {QUAD_GAUSS,QUAD_LOBATTO} QuadratureType;
87: static const char *QuadratureTypes[] = {"gauss","lobatto","QuadratureType","QUAD_",0};
88: PETSC_UNUSED static const PetscReal HexQWeights[8] = {1,1,1,1,1,1,1,1};
89: PETSC_UNUSED static const PetscReal HexQNodes[] = {-0.57735026918962573, 0.57735026918962573};
90: #define G 0.57735026918962573
91: #define H (0.5*(1.+G))
92: #define L (0.5*(1.-G))
93: #define M (-0.5)
94: #define P (0.5)
95: /* Special quadrature: Lobatto in horizontal, Gauss in vertical */
96: static const PetscReal HexQInterp_Lobatto[8][8] = {{H,0,0,0,L,0,0,0},
97: {0,H,0,0,0,L,0,0},
98: {0,0,H,0,0,0,L,0},
99: {0,0,0,H,0,0,0,L},
100: {L,0,0,0,H,0,0,0},
101: {0,L,0,0,0,H,0,0},
102: {0,0,L,0,0,0,H,0},
103: {0,0,0,L,0,0,0,H}};
104: static const PetscReal HexQDeriv_Lobatto[8][8][3] = {
105: {{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} },
106: {{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} },
107: {{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} },
108: {{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}},
109: {{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} },
110: {{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} },
111: {{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} },
112: {{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}}};
113: /* Stanndard Gauss */
114: 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},
115: {L*H*H,H*H*H,H*L*H,L*L*H, L*H*L,H*H*L,H*L*L,L*L*L},
116: {L*L*H,H*L*H,H*H*H,L*H*H, L*L*L,H*L*L,H*H*L,L*H*L},
117: {H*L*H,L*L*H,L*H*H,H*H*H, H*L*L,L*L*L,L*H*L,H*H*L},
118: {H*H*L,L*H*L,L*L*L,H*L*L, H*H*H,L*H*H,L*L*H,H*L*H},
119: {L*H*L,H*H*L,H*L*L,L*L*L, L*H*H,H*H*H,H*L*H,L*L*H},
120: {L*L*L,H*L*L,H*H*L,L*H*L, L*L*H,H*L*H,H*H*H,L*H*H},
121: {H*L*L,L*L*L,L*H*L,H*H*L, H*L*H,L*L*H,L*H*H,H*H*H}};
122: static const PetscReal HexQDeriv_Gauss[8][8][3] = {
123: {{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}},
124: {{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}},
125: {{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}},
126: {{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}},
127: {{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}},
128: {{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}},
129: {{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}},
130: {{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}}};
131: static const PetscReal (*HexQInterp)[8],(*HexQDeriv)[8][3];
132: /* Standard 2x2 Gauss quadrature for the bottom layer. */
133: static const PetscReal QuadQInterp[4][4] = {{H*H,L*H,L*L,H*L},
134: {L*H,H*H,H*L,L*L},
135: {L*L,H*L,H*H,L*H},
136: {H*L,L*L,L*H,H*H}};
137: static const PetscReal QuadQDeriv[4][4][2] = {
138: {{M*H,M*H},{P*H,M*L},{P*L,P*L},{M*L,P*H}},
139: {{M*H,M*L},{P*H,M*H},{P*L,P*H},{M*L,P*L}},
140: {{M*L,M*L},{P*L,M*H},{P*H,P*H},{M*H,P*L}},
141: {{M*L,M*H},{P*L,M*L},{P*H,P*L},{M*H,P*H}}};
142: #undef G
143: #undef H
144: #undef L
145: #undef M
146: #undef P
148: #define HexExtract(x,i,j,k,n) do { \
149: (n)[0] = (x)[i][j][k]; \
150: (n)[1] = (x)[i+1][j][k]; \
151: (n)[2] = (x)[i+1][j+1][k]; \
152: (n)[3] = (x)[i][j+1][k]; \
153: (n)[4] = (x)[i][j][k+1]; \
154: (n)[5] = (x)[i+1][j][k+1]; \
155: (n)[6] = (x)[i+1][j+1][k+1]; \
156: (n)[7] = (x)[i][j+1][k+1]; \
157: } while (0)
159: #define HexExtractRef(x,i,j,k,n) do { \
160: (n)[0] = &(x)[i][j][k]; \
161: (n)[1] = &(x)[i+1][j][k]; \
162: (n)[2] = &(x)[i+1][j+1][k]; \
163: (n)[3] = &(x)[i][j+1][k]; \
164: (n)[4] = &(x)[i][j][k+1]; \
165: (n)[5] = &(x)[i+1][j][k+1]; \
166: (n)[6] = &(x)[i+1][j+1][k+1]; \
167: (n)[7] = &(x)[i][j+1][k+1]; \
168: } while (0)
170: #define QuadExtract(x,i,j,n) do { \
171: (n)[0] = (x)[i][j]; \
172: (n)[1] = (x)[i+1][j]; \
173: (n)[2] = (x)[i+1][j+1]; \
174: (n)[3] = (x)[i][j+1]; \
175: } while (0)
177: static void HexGrad(const PetscReal dphi[][3],const PetscReal zn[],PetscReal dz[])
178: {
179: PetscInt i;
180: dz[0] = dz[1] = dz[2] = 0;
181: for (i=0; i<8; i++) {
182: dz[0] += dphi[i][0] * zn[i];
183: dz[1] += dphi[i][1] * zn[i];
184: dz[2] += dphi[i][2] * zn[i];
185: }
186: }
188: static void HexComputeGeometry(PetscInt q,PetscReal hx,PetscReal hy,const PetscReal dz[PETSC_RESTRICT],PetscReal phi[PETSC_RESTRICT],PetscReal dphi[PETSC_RESTRICT][3],PetscReal *PETSC_RESTRICT jw)
189: {
190: const PetscReal jac[3][3] = {{hx/2,0,0}, {0,hy/2,0}, {dz[0],dz[1],dz[2]}};
191: const PetscReal 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]}};
192: const PetscReal jdet = jac[0][0]*jac[1][1]*jac[2][2];
193: PetscInt i;
195: for (i=0; i<8; i++) {
196: const PetscReal *dphir = HexQDeriv[q][i];
197: phi[i] = HexQInterp[q][i];
198: dphi[i][0] = dphir[0]*ijac[0][0] + dphir[1]*ijac[1][0] + dphir[2]*ijac[2][0];
199: dphi[i][1] = dphir[0]*ijac[0][1] + dphir[1]*ijac[1][1] + dphir[2]*ijac[2][1];
200: dphi[i][2] = dphir[0]*ijac[0][2] + dphir[1]*ijac[1][2] + dphir[2]*ijac[2][2];
201: }
202: *jw = 1.0 * jdet;
203: }
205: typedef struct _p_THI *THI;
206: typedef struct _n_Units *Units;
208: typedef struct {
209: PetscScalar u,v;
210: } Node;
212: typedef struct {
213: PetscScalar b; /* bed */
214: PetscScalar h; /* thickness */
215: PetscScalar beta2; /* friction */
216: } PrmNode;
218: typedef struct {
219: PetscReal min,max,cmin,cmax;
220: } PRange;
222: typedef enum {THIASSEMBLY_TRIDIAGONAL,THIASSEMBLY_FULL} THIAssemblyMode;
224: struct _p_THI {
225: PETSCHEADER(int);
226: void (*initialize)(THI,PetscReal x,PetscReal y,PrmNode *p);
227: PetscInt zlevels;
228: PetscReal Lx,Ly,Lz; /* Model domain */
229: PetscReal alpha; /* Bed angle */
230: Units units;
231: PetscReal dirichlet_scale;
232: PetscReal ssa_friction_scale;
233: PRange eta;
234: PRange beta2;
235: struct {
236: PetscReal Bd2,eps,exponent;
237: } viscosity;
238: struct {
239: PetscReal irefgam,eps2,exponent,refvel,epsvel;
240: } friction;
241: PetscReal rhog;
242: PetscBool no_slip;
243: PetscBool tridiagonal;
244: PetscBool coarse2d;
245: PetscBool verbose;
246: MatType mattype;
247: };
249: struct _n_Units {
250: /* fundamental */
251: PetscReal meter;
252: PetscReal kilogram;
253: PetscReal second;
254: /* derived */
255: PetscReal Pascal;
256: PetscReal year;
257: };
259: static PetscErrorCode THIJacobianLocal_3D_Full(DMDALocalInfo*,Node***,Mat,Mat,THI);
260: static PetscErrorCode THIJacobianLocal_3D_Tridiagonal(DMDALocalInfo*,Node***,Mat,Mat,THI);
261: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo*,Node**,Mat,Mat,THI);
263: static void PrmHexGetZ(const PrmNode pn[],PetscInt k,PetscInt zm,PetscReal zn[])
264: {
265: const PetscScalar zm1 = zm-1,
266: znl[8] = {pn[0].b + pn[0].h*(PetscScalar)k/zm1,
267: pn[1].b + pn[1].h*(PetscScalar)k/zm1,
268: pn[2].b + pn[2].h*(PetscScalar)k/zm1,
269: pn[3].b + pn[3].h*(PetscScalar)k/zm1,
270: pn[0].b + pn[0].h*(PetscScalar)(k+1)/zm1,
271: pn[1].b + pn[1].h*(PetscScalar)(k+1)/zm1,
272: pn[2].b + pn[2].h*(PetscScalar)(k+1)/zm1,
273: pn[3].b + pn[3].h*(PetscScalar)(k+1)/zm1};
274: PetscInt i;
275: for (i=0; i<8; i++) zn[i] = PetscRealPart(znl[i]);
276: }
278: /* Tests A and C are from the ISMIP-HOM paper (Pattyn et al. 2008) */
279: static void THIInitialize_HOM_A(THI thi,PetscReal x,PetscReal y,PrmNode *p)
280: {
281: Units units = thi->units;
282: PetscReal s = -x*PetscSinReal(thi->alpha);
284: p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx) * PetscSinReal(y*2*PETSC_PI/thi->Ly);
285: p->h = s - p->b;
286: p->beta2 = 1e30;
287: }
289: static void THIInitialize_HOM_C(THI thi,PetscReal x,PetscReal y,PrmNode *p)
290: {
291: Units units = thi->units;
292: PetscReal s = -x*PetscSinReal(thi->alpha);
294: p->b = s - 1000*units->meter;
295: p->h = s - p->b;
296: /* tau_b = beta2 v is a stress (Pa) */
297: p->beta2 = 1000 * (1 + PetscSinReal(x*2*PETSC_PI/thi->Lx)*PetscSinReal(y*2*PETSC_PI/thi->Ly)) * units->Pascal * units->year / units->meter;
298: }
300: /* These are just toys */
302: /* Same bed as test A, free slip everywhere except for a discontinuous jump to a circular sticky region in the middle. */
303: static void THIInitialize_HOM_X(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
304: {
305: Units units = thi->units;
306: PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
307: PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
308: p->b = s - 1000*units->meter + 500*units->meter*PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
309: p->h = s - p->b;
310: p->beta2 = 1000 * (r < 1 ? 2 : 0) * units->Pascal * units->year / units->meter;
311: }
313: /* Like Z, but with 200 meter cliffs */
314: static void THIInitialize_HOM_Y(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
315: {
316: Units units = thi->units;
317: PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
318: PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
320: p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
321: if (PetscRealPart(p->b) > -700*units->meter) p->b += 200*units->meter;
322: p->h = s - p->b;
323: 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;
324: }
326: /* Same bed as A, smoothly varying slipperiness, similar to MATLAB's "sombrero" (uncorrelated with bathymetry) */
327: static void THIInitialize_HOM_Z(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
328: {
329: Units units = thi->units;
330: PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
331: PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
333: p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
334: p->h = s - p->b;
335: 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;
336: }
338: static void THIFriction(THI thi,PetscReal rbeta2,PetscReal gam,PetscReal *beta2,PetscReal *dbeta2)
339: {
340: if (thi->friction.irefgam == 0) {
341: Units units = thi->units;
342: thi->friction.irefgam = 1./(0.5*PetscSqr(thi->friction.refvel * units->meter / units->year));
343: thi->friction.eps2 = 0.5*PetscSqr(thi->friction.epsvel * units->meter / units->year) * thi->friction.irefgam;
344: }
345: if (thi->friction.exponent == 0) {
346: *beta2 = rbeta2;
347: *dbeta2 = 0;
348: } else {
349: *beta2 = rbeta2 * PetscPowReal(thi->friction.eps2 + gam*thi->friction.irefgam,thi->friction.exponent);
350: *dbeta2 = thi->friction.exponent * *beta2 / (thi->friction.eps2 + gam*thi->friction.irefgam) * thi->friction.irefgam;
351: }
352: }
354: static void THIViscosity(THI thi,PetscReal gam,PetscReal *eta,PetscReal *deta)
355: {
356: PetscReal Bd2,eps,exponent;
357: if (thi->viscosity.Bd2 == 0) {
358: Units units = thi->units;
359: const PetscReal
360: n = 3., /* Glen exponent */
361: p = 1. + 1./n, /* for Stokes */
362: A = 1.e-16 * PetscPowReal(units->Pascal,-n) / units->year, /* softness parameter (Pa^{-n}/s) */
363: B = PetscPowReal(A,-1./n); /* hardness parameter */
364: thi->viscosity.Bd2 = B/2;
365: thi->viscosity.exponent = (p-2)/2;
366: thi->viscosity.eps = 0.5*PetscSqr(1e-5 / units->year);
367: }
368: Bd2 = thi->viscosity.Bd2;
369: exponent = thi->viscosity.exponent;
370: eps = thi->viscosity.eps;
371: *eta = Bd2 * PetscPowReal(eps + gam,exponent);
372: *deta = exponent * (*eta) / (eps + gam);
373: }
375: static void RangeUpdate(PetscReal *min,PetscReal *max,PetscReal x)
376: {
377: if (x < *min) *min = x;
378: if (x > *max) *max = x;
379: }
381: static void PRangeClear(PRange *p)
382: {
383: p->cmin = p->min = 1e100;
384: p->cmax = p->max = -1e100;
385: }
389: static PetscErrorCode PRangeMinMax(PRange *p,PetscReal min,PetscReal max)
390: {
393: p->cmin = min;
394: p->cmax = max;
395: if (min < p->min) p->min = min;
396: if (max > p->max) p->max = max;
397: return(0);
398: }
402: static PetscErrorCode THIDestroy(THI *thi)
403: {
407: if (!*thi) return(0);
408: if (--((PetscObject)(*thi))->refct > 0) {*thi = 0; return(0);}
409: PetscFree((*thi)->units);
410: PetscFree((*thi)->mattype);
411: PetscHeaderDestroy(thi);
412: return(0);
413: }
417: static PetscErrorCode THICreate(MPI_Comm comm,THI *inthi)
418: {
419: static PetscBool registered = PETSC_FALSE;
420: THI thi;
421: Units units;
422: PetscErrorCode ierr;
425: *inthi = 0;
426: if (!registered) {
427: PetscClassIdRegister("Toy Hydrostatic Ice",&THI_CLASSID);
428: registered = PETSC_TRUE;
429: }
430: PetscHeaderCreate(thi,THI_CLASSID,"THI","Toy Hydrostatic Ice","",comm,THIDestroy,0);
432: PetscNew(&thi->units);
433: units = thi->units;
434: units->meter = 1e-2;
435: units->second = 1e-7;
436: units->kilogram = 1e-12;
438: PetscOptionsBegin(comm,NULL,"Scaled units options","");
439: {
440: PetscOptionsReal("-units_meter","1 meter in scaled length units","",units->meter,&units->meter,NULL);
441: PetscOptionsReal("-units_second","1 second in scaled time units","",units->second,&units->second,NULL);
442: PetscOptionsReal("-units_kilogram","1 kilogram in scaled mass units","",units->kilogram,&units->kilogram,NULL);
443: }
444: PetscOptionsEnd();
445: units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second));
446: units->year = 31556926. * units->second, /* seconds per year */
448: thi->Lx = 10.e3;
449: thi->Ly = 10.e3;
450: thi->Lz = 1000;
451: thi->dirichlet_scale = 1;
452: thi->verbose = PETSC_FALSE;
454: PetscOptionsBegin(comm,NULL,"Toy Hydrostatic Ice options","");
455: {
456: QuadratureType quad = QUAD_GAUSS;
457: char homexp[] = "A";
458: char mtype[256] = MATSBAIJ;
459: PetscReal L,m = 1.0;
460: PetscBool flg;
461: L = thi->Lx;
462: PetscOptionsReal("-thi_L","Domain size (m)","",L,&L,&flg);
463: if (flg) thi->Lx = thi->Ly = L;
464: PetscOptionsReal("-thi_Lx","X Domain size (m)","",thi->Lx,&thi->Lx,NULL);
465: PetscOptionsReal("-thi_Ly","Y Domain size (m)","",thi->Ly,&thi->Ly,NULL);
466: PetscOptionsReal("-thi_Lz","Z Domain size (m)","",thi->Lz,&thi->Lz,NULL);
467: PetscOptionsString("-thi_hom","ISMIP-HOM experiment (A or C)","",homexp,homexp,sizeof(homexp),NULL);
468: switch (homexp[0] = toupper(homexp[0])) {
469: case 'A':
470: thi->initialize = THIInitialize_HOM_A;
471: thi->no_slip = PETSC_TRUE;
472: thi->alpha = 0.5;
473: break;
474: case 'C':
475: thi->initialize = THIInitialize_HOM_C;
476: thi->no_slip = PETSC_FALSE;
477: thi->alpha = 0.1;
478: break;
479: case 'X':
480: thi->initialize = THIInitialize_HOM_X;
481: thi->no_slip = PETSC_FALSE;
482: thi->alpha = 0.3;
483: break;
484: case 'Y':
485: thi->initialize = THIInitialize_HOM_Y;
486: thi->no_slip = PETSC_FALSE;
487: thi->alpha = 0.5;
488: break;
489: case 'Z':
490: thi->initialize = THIInitialize_HOM_Z;
491: thi->no_slip = PETSC_FALSE;
492: thi->alpha = 0.5;
493: break;
494: default:
495: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"HOM experiment '%c' not implemented",homexp[0]);
496: }
497: PetscOptionsEnum("-thi_quadrature","Quadrature to use for 3D elements","",QuadratureTypes,(PetscEnum)quad,(PetscEnum*)&quad,NULL);
498: switch (quad) {
499: case QUAD_GAUSS:
500: HexQInterp = HexQInterp_Gauss;
501: HexQDeriv = HexQDeriv_Gauss;
502: break;
503: case QUAD_LOBATTO:
504: HexQInterp = HexQInterp_Lobatto;
505: HexQDeriv = HexQDeriv_Lobatto;
506: break;
507: }
508: PetscOptionsReal("-thi_alpha","Bed angle (degrees)","",thi->alpha,&thi->alpha,NULL);
510: thi->friction.refvel = 100.;
511: thi->friction.epsvel = 1.;
513: PetscOptionsReal("-thi_friction_refvel","Reference velocity for sliding","",thi->friction.refvel,&thi->friction.refvel,NULL);
514: PetscOptionsReal("-thi_friction_epsvel","Regularization velocity for sliding","",thi->friction.epsvel,&thi->friction.epsvel,NULL);
515: PetscOptionsReal("-thi_friction_m","Friction exponent, 0=Coulomb, 1=Navier","",m,&m,NULL);
517: thi->friction.exponent = (m-1)/2;
519: PetscOptionsReal("-thi_dirichlet_scale","Scale Dirichlet boundary conditions by this factor","",thi->dirichlet_scale,&thi->dirichlet_scale,NULL);
520: 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);
521: PetscOptionsBool("-thi_coarse2d","Use a 2D coarse space corresponding to SSA","",thi->coarse2d,&thi->coarse2d,NULL);
522: PetscOptionsBool("-thi_tridiagonal","Assemble a tridiagonal system (column coupling only) on the finest level","",thi->tridiagonal,&thi->tridiagonal,NULL);
523: PetscOptionsFList("-thi_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)mtype,sizeof(mtype),NULL);
524: PetscStrallocpy(mtype,(char**)&thi->mattype);
525: PetscOptionsBool("-thi_verbose","Enable verbose output (like matrix sizes and statistics)","",thi->verbose,&thi->verbose,NULL);
526: }
527: PetscOptionsEnd();
529: /* dimensionalize */
530: thi->Lx *= units->meter;
531: thi->Ly *= units->meter;
532: thi->Lz *= units->meter;
533: thi->alpha *= PETSC_PI / 180;
535: PRangeClear(&thi->eta);
536: PRangeClear(&thi->beta2);
538: {
539: PetscReal u = 1000*units->meter/(3e7*units->second),
540: gradu = u / (100*units->meter),eta,deta,
541: rho = 910 * units->kilogram/PetscPowReal(units->meter,3),
542: grav = 9.81 * units->meter/PetscSqr(units->second),
543: driving = rho * grav * PetscSinReal(thi->alpha) * 1000*units->meter;
544: THIViscosity(thi,0.5*gradu*gradu,&eta,&deta);
545: thi->rhog = rho * grav;
546: if (thi->verbose) {
547: PetscPrintf(PetscObjectComm((PetscObject)thi),"Units: meter %8.2g second %8.2g kg %8.2g Pa %8.2g\n",(double)units->meter,(double)units->second,(double)units->kilogram,(double)units->Pascal);
548: PetscPrintf(PetscObjectComm((PetscObject)thi),"Domain (%6.2g,%6.2g,%6.2g), pressure %8.2g, driving stress %8.2g\n",(double)thi->Lx,(double)thi->Ly,(double)thi->Lz,(double)(rho*grav*1e3*units->meter),(double)driving);
549: PetscPrintf(PetscObjectComm((PetscObject)thi),"Large velocity 1km/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",(double)u,(double)gradu,(double)eta,(double)(2*eta*gradu),(double)(2*eta*gradu/driving));
550: THIViscosity(thi,0.5*PetscSqr(1e-3*gradu),&eta,&deta);
551: PetscPrintf(PetscObjectComm((PetscObject)thi),"Small velocity 1m/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",(double)(1e-3*u),(double)(1e-3*gradu),(double)eta,(double)(2*eta*1e-3*gradu),(double)(2*eta*1e-3*gradu/driving));
552: }
553: }
555: *inthi = thi;
556: return(0);
557: }
561: static PetscErrorCode THIInitializePrm(THI thi,DM da2prm,Vec prm)
562: {
563: PrmNode **p;
564: PetscInt i,j,xs,xm,ys,ym,mx,my;
568: DMDAGetGhostCorners(da2prm,&ys,&xs,0,&ym,&xm,0);
569: DMDAGetInfo(da2prm,0, &my,&mx,0, 0,0,0, 0,0,0,0,0,0);
570: DMDAVecGetArray(da2prm,prm,&p);
571: for (i=xs; i<xs+xm; i++) {
572: for (j=ys; j<ys+ym; j++) {
573: PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my;
574: thi->initialize(thi,xx,yy,&p[i][j]);
575: }
576: }
577: DMDAVecRestoreArray(da2prm,prm,&p);
578: return(0);
579: }
583: static PetscErrorCode THISetUpDM(THI thi,DM dm)
584: {
585: PetscErrorCode ierr;
586: PetscInt refinelevel,coarsenlevel,level,dim,Mx,My,Mz,mx,my,s;
587: DMDAStencilType st;
588: DM da2prm;
589: Vec X;
592: DMDAGetInfo(dm,&dim, &Mz,&My,&Mx, 0,&my,&mx, 0,&s,0,0,0,&st);
593: if (dim == 2) {
594: DMDAGetInfo(dm,&dim, &My,&Mx,0, &my,&mx,0, 0,&s,0,0,0,&st);
595: }
596: DMGetRefineLevel(dm,&refinelevel);
597: DMGetCoarsenLevel(dm,&coarsenlevel);
598: level = refinelevel - coarsenlevel;
599: DMDACreate2d(PetscObjectComm((PetscObject)thi),DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,My,Mx,my,mx,sizeof(PrmNode)/sizeof(PetscScalar),s,0,0,&da2prm);
600: DMCreateLocalVector(da2prm,&X);
601: {
602: PetscReal Lx = thi->Lx / thi->units->meter,Ly = thi->Ly / thi->units->meter,Lz = thi->Lz / thi->units->meter;
603: if (dim == 2) {
604: PetscPrintf(PetscObjectComm((PetscObject)thi),"Level %D domain size (m) %8.2g x %8.2g, num elements %D x %D (%D), size (m) %g x %g\n",level,(double)Lx,(double)Ly,Mx,My,Mx*My,(double)(Lx/Mx),(double)(Ly/My));
605: } else {
606: PetscPrintf(PetscObjectComm((PetscObject)thi),"Level %D domain size (m) %8.2g x %8.2g x %8.2g, num elements %D x %D x %D (%D), size (m) %g x %g x %g\n",level,(double)Lx,(double)Ly,(double)Lz,Mx,My,Mz,Mx*My*Mz,(double)(Lx/Mx),(double)(Ly/My),(double)(1000./(Mz-1)));
607: }
608: }
609: THIInitializePrm(thi,da2prm,X);
610: if (thi->tridiagonal) { /* Reset coarse Jacobian evaluation */
611: DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_3D_Full,thi);
612: }
613: if (thi->coarse2d) {
614: DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_2D,thi);
615: }
616: PetscObjectCompose((PetscObject)dm,"DMDA2Prm",(PetscObject)da2prm);
617: PetscObjectCompose((PetscObject)dm,"DMDA2Prm_Vec",(PetscObject)X);
618: DMDestroy(&da2prm);
619: VecDestroy(&X);
620: return(0);
621: }
625: static PetscErrorCode DMCoarsenHook_THI(DM dmf,DM dmc,void *ctx)
626: {
627: THI thi = (THI)ctx;
629: PetscInt rlevel,clevel;
632: THISetUpDM(thi,dmc);
633: DMGetRefineLevel(dmc,&rlevel);
634: DMGetCoarsenLevel(dmc,&clevel);
635: if (rlevel-clevel == 0) {DMSetMatType(dmc,MATAIJ);}
636: DMCoarsenHookAdd(dmc,DMCoarsenHook_THI,NULL,thi);
637: return(0);
638: }
642: static PetscErrorCode DMRefineHook_THI(DM dmc,DM dmf,void *ctx)
643: {
644: THI thi = (THI)ctx;
648: THISetUpDM(thi,dmf);
649: DMSetMatType(dmf,thi->mattype);
650: DMRefineHookAdd(dmf,DMRefineHook_THI,NULL,thi);
651: /* With grid sequencing, a formerly-refined DM will later be coarsened by PCSetUp_MG */
652: DMCoarsenHookAdd(dmf,DMCoarsenHook_THI,NULL,thi);
653: return(0);
654: }
658: static PetscErrorCode THIDAGetPrm(DM da,PrmNode ***prm)
659: {
661: DM da2prm;
662: Vec X;
665: PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
666: if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
667: PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
668: if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
669: DMDAVecGetArray(da2prm,X,prm);
670: return(0);
671: }
675: static PetscErrorCode THIDARestorePrm(DM da,PrmNode ***prm)
676: {
678: DM da2prm;
679: Vec X;
682: PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
683: if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
684: PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
685: if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
686: DMDAVecRestoreArray(da2prm,X,prm);
687: return(0);
688: }
692: static PetscErrorCode THIInitial(SNES snes,Vec X,void *ctx)
693: {
694: THI thi;
695: PetscInt i,j,k,xs,xm,ys,ym,zs,zm,mx,my;
696: PetscReal hx,hy;
697: PrmNode **prm;
698: Node ***x;
700: DM da;
703: SNESGetDM(snes,&da);
704: DMGetApplicationContext(da,&thi);
705: DMDAGetInfo(da,0, 0,&my,&mx, 0,0,0, 0,0,0,0,0,0);
706: DMDAGetCorners(da,&zs,&ys,&xs,&zm,&ym,&xm);
707: DMDAVecGetArray(da,X,&x);
708: THIDAGetPrm(da,&prm);
709: hx = thi->Lx / mx;
710: hy = thi->Ly / my;
711: for (i=xs; i<xs+xm; i++) {
712: for (j=ys; j<ys+ym; j++) {
713: for (k=zs; k<zs+zm; k++) {
714: const PetscScalar zm1 = zm-1,
715: 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),
716: 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);
717: x[i][j][k].u = 0. * drivingx * prm[i][j].h*(PetscScalar)k/zm1;
718: x[i][j][k].v = 0. * drivingy * prm[i][j].h*(PetscScalar)k/zm1;
719: }
720: }
721: }
722: DMDAVecRestoreArray(da,X,&x);
723: THIDARestorePrm(da,&prm);
724: return(0);
725: }
727: static void PointwiseNonlinearity(THI thi,const Node n[PETSC_RESTRICT],const PetscReal phi[PETSC_RESTRICT],PetscReal dphi[PETSC_RESTRICT][3],PetscScalar *PETSC_RESTRICT u,PetscScalar *PETSC_RESTRICT v,PetscScalar du[PETSC_RESTRICT],PetscScalar dv[PETSC_RESTRICT],PetscReal *eta,PetscReal *deta)
728: {
729: PetscInt l,ll;
730: PetscScalar gam;
732: du[0] = du[1] = du[2] = 0;
733: dv[0] = dv[1] = dv[2] = 0;
734: *u = 0;
735: *v = 0;
736: for (l=0; l<8; l++) {
737: *u += phi[l] * n[l].u;
738: *v += phi[l] * n[l].v;
739: for (ll=0; ll<3; ll++) {
740: du[ll] += dphi[l][ll] * n[l].u;
741: dv[ll] += dphi[l][ll] * n[l].v;
742: }
743: }
744: gam = PetscSqr(du[0]) + PetscSqr(dv[1]) + du[0]*dv[1] + 0.25*PetscSqr(du[1]+dv[0]) + 0.25*PetscSqr(du[2]) + 0.25*PetscSqr(dv[2]);
745: THIViscosity(thi,PetscRealPart(gam),eta,deta);
746: }
748: static void PointwiseNonlinearity2D(THI thi,Node n[],PetscReal phi[],PetscReal dphi[4][2],PetscScalar *u,PetscScalar *v,PetscScalar du[],PetscScalar dv[],PetscReal *eta,PetscReal *deta)
749: {
750: PetscInt l,ll;
751: PetscScalar gam;
753: du[0] = du[1] = 0;
754: dv[0] = dv[1] = 0;
755: *u = 0;
756: *v = 0;
757: for (l=0; l<4; l++) {
758: *u += phi[l] * n[l].u;
759: *v += phi[l] * n[l].v;
760: for (ll=0; ll<2; ll++) {
761: du[ll] += dphi[l][ll] * n[l].u;
762: dv[ll] += dphi[l][ll] * n[l].v;
763: }
764: }
765: gam = PetscSqr(du[0]) + PetscSqr(dv[1]) + du[0]*dv[1] + 0.25*PetscSqr(du[1]+dv[0]);
766: THIViscosity(thi,PetscRealPart(gam),eta,deta);
767: }
771: static PetscErrorCode THIFunctionLocal(DMDALocalInfo *info,Node ***x,Node ***f,THI thi)
772: {
773: PetscInt xs,ys,xm,ym,zm,i,j,k,q,l;
774: PetscReal hx,hy,etamin,etamax,beta2min,beta2max;
775: PrmNode **prm;
779: xs = info->zs;
780: ys = info->ys;
781: xm = info->zm;
782: ym = info->ym;
783: zm = info->xm;
784: hx = thi->Lx / info->mz;
785: hy = thi->Ly / info->my;
787: etamin = 1e100;
788: etamax = 0;
789: beta2min = 1e100;
790: beta2max = 0;
792: THIDAGetPrm(info->da,&prm);
794: for (i=xs; i<xs+xm; i++) {
795: for (j=ys; j<ys+ym; j++) {
796: PrmNode pn[4];
797: QuadExtract(prm,i,j,pn);
798: for (k=0; k<zm-1; k++) {
799: PetscInt ls = 0;
800: Node n[8],*fn[8];
801: PetscReal zn[8],etabase = 0;
802: PrmHexGetZ(pn,k,zm,zn);
803: HexExtract(x,i,j,k,n);
804: HexExtractRef(f,i,j,k,fn);
805: if (thi->no_slip && k == 0) {
806: for (l=0; l<4; l++) n[l].u = n[l].v = 0;
807: /* The first 4 basis functions lie on the bottom layer, so their contribution is exactly 0, hence we can skip them */
808: ls = 4;
809: }
810: for (q=0; q<8; q++) {
811: PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta;
812: PetscScalar du[3],dv[3],u,v;
813: HexGrad(HexQDeriv[q],zn,dz);
814: HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
815: PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
816: jw /= thi->rhog; /* scales residuals to be O(1) */
817: if (q == 0) etabase = eta;
818: RangeUpdate(&etamin,&etamax,eta);
819: for (l=ls; l<8; l++) { /* test functions */
820: const PetscReal ds[2] = {-PetscSinReal(thi->alpha),0};
821: const PetscReal pp = phi[l],*dp = dphi[l];
822: 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];
823: 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];
824: }
825: }
826: if (k == 0) { /* we are on a bottom face */
827: if (thi->no_slip) {
828: /* Note: Non-Galerkin coarse grid operators are very sensitive to the scaling of Dirichlet boundary
829: * conditions. After shenanigans above, etabase contains the effective viscosity at the closest quadrature
830: * point to the bed. We want the diagonal entry in the Dirichlet condition to have similar magnitude to the
831: * diagonal entry corresponding to the adjacent node. The fundamental scaling of the viscous part is in
832: * diagu, diagv below. This scaling is easy to recognize by considering the finite difference operator after
833: * scaling by element size. The no-slip Dirichlet condition is scaled by this factor, and also in the
834: * assembled matrix (see the similar block in THIJacobianLocal).
835: *
836: * Note that the residual at this Dirichlet node is linear in the state at this node, but also depends
837: * (nonlinearly in general) on the neighboring interior nodes through the local viscosity. This will make
838: * a matrix-free Jacobian have extra entries in the corresponding row. We assemble only the diagonal part,
839: * so the solution will exactly satisfy the boundary condition after the first linear iteration.
840: */
841: const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1.);
842: 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);
843: fn[0]->u = thi->dirichlet_scale*diagu*x[i][j][k].u;
844: fn[0]->v = thi->dirichlet_scale*diagv*x[i][j][k].v;
845: } else { /* Integrate over bottom face to apply boundary condition */
846: for (q=0; q<4; q++) {
847: const PetscReal jw = 0.25*hx*hy/thi->rhog,*phi = QuadQInterp[q];
848: PetscScalar u =0,v=0,rbeta2=0;
849: PetscReal beta2,dbeta2;
850: for (l=0; l<4; l++) {
851: u += phi[l]*n[l].u;
852: v += phi[l]*n[l].v;
853: rbeta2 += phi[l]*pn[l].beta2;
854: }
855: THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
856: RangeUpdate(&beta2min,&beta2max,beta2);
857: for (l=0; l<4; l++) {
858: const PetscReal pp = phi[l];
859: fn[ls+l]->u += pp*jw*beta2*u;
860: fn[ls+l]->v += pp*jw*beta2*v;
861: }
862: }
863: }
864: }
865: }
866: }
867: }
869: THIDARestorePrm(info->da,&prm);
871: PRangeMinMax(&thi->eta,etamin,etamax);
872: PRangeMinMax(&thi->beta2,beta2min,beta2max);
873: return(0);
874: }
878: static PetscErrorCode THIMatrixStatistics(THI thi,Mat B,PetscViewer viewer)
879: {
881: PetscReal nrm;
882: PetscInt m;
883: PetscMPIInt rank;
886: MatNorm(B,NORM_FROBENIUS,&nrm);
887: MatGetSize(B,&m,0);
888: MPI_Comm_rank(PetscObjectComm((PetscObject)B),&rank);
889: if (!rank) {
890: PetscScalar val0,val2;
891: MatGetValue(B,0,0,&val0);
892: MatGetValue(B,2,2,&val2);
893: PetscViewerASCIIPrintf(viewer,"Matrix dim %D norm %8.2e (0,0) %8.2e (2,2) %8.2e %8.2e <= eta <= %8.2e %8.2e <= beta2 <= %8.2e\n",m,(double)nrm,(double)PetscRealPart(val0),(double)PetscRealPart(val2),(double)thi->eta.cmin,(double)thi->eta.cmax,(double)thi->beta2.cmin,(double)thi->beta2.cmax);
894: }
895: return(0);
896: }
900: static PetscErrorCode THISurfaceStatistics(DM da,Vec X,PetscReal *min,PetscReal *max,PetscReal *mean)
901: {
903: Node ***x;
904: PetscInt i,j,xs,ys,zs,xm,ym,zm,mx,my,mz;
905: PetscReal umin = 1e100,umax=-1e100;
906: PetscScalar usum = 0.0,gusum;
909: *min = *max = *mean = 0;
910: DMDAGetInfo(da,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
911: DMDAGetCorners(da,&zs,&ys,&xs,&zm,&ym,&xm);
912: if (zs != 0 || zm != mz) SETERRQ(PETSC_COMM_SELF,1,"Unexpected decomposition");
913: DMDAVecGetArray(da,X,&x);
914: for (i=xs; i<xs+xm; i++) {
915: for (j=ys; j<ys+ym; j++) {
916: PetscReal u = PetscRealPart(x[i][j][zm-1].u);
917: RangeUpdate(&umin,&umax,u);
918: usum += u;
919: }
920: }
921: DMDAVecRestoreArray(da,X,&x);
922: MPI_Allreduce(&umin,min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));
923: MPI_Allreduce(&umax,max,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));
924: MPI_Allreduce(&usum,&gusum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da));
925: *mean = PetscRealPart(gusum) / (mx*my);
926: return(0);
927: }
931: static PetscErrorCode THISolveStatistics(THI thi,SNES snes,PetscInt coarsened,const char name[])
932: {
933: MPI_Comm comm;
934: Vec X;
935: DM dm;
939: PetscObjectGetComm((PetscObject)thi,&comm);
940: SNESGetSolution(snes,&X);
941: SNESGetDM(snes,&dm);
942: PetscPrintf(comm,"Solution statistics after solve: %s\n",name);
943: {
944: PetscInt its,lits;
945: SNESConvergedReason reason;
946: SNESGetIterationNumber(snes,&its);
947: SNESGetConvergedReason(snes,&reason);
948: SNESGetLinearSolveIterations(snes,&lits);
949: PetscPrintf(comm,"%s: Number of SNES iterations = %D, total linear iterations = %D\n",SNESConvergedReasons[reason],its,lits);
950: }
951: {
952: PetscReal nrm2,tmin[3]={1e100,1e100,1e100},tmax[3]={-1e100,-1e100,-1e100},min[3],max[3];
953: PetscInt i,j,m;
954: const PetscScalar *x;
955: VecNorm(X,NORM_2,&nrm2);
956: VecGetLocalSize(X,&m);
957: VecGetArrayRead(X,&x);
958: for (i=0; i<m; i+=2) {
959: PetscReal u = PetscRealPart(x[i]),v = PetscRealPart(x[i+1]),c = PetscSqrtReal(u*u+v*v);
960: tmin[0] = PetscMin(u,tmin[0]);
961: tmin[1] = PetscMin(v,tmin[1]);
962: tmin[2] = PetscMin(c,tmin[2]);
963: tmax[0] = PetscMax(u,tmax[0]);
964: tmax[1] = PetscMax(v,tmax[1]);
965: tmax[2] = PetscMax(c,tmax[2]);
966: }
967: VecRestoreArrayRead(X,&x);
968: MPI_Allreduce(tmin,min,3,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)thi));
969: MPI_Allreduce(tmax,max,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)thi));
970: /* Dimensionalize to meters/year */
971: nrm2 *= thi->units->year / thi->units->meter;
972: for (j=0; j<3; j++) {
973: min[j] *= thi->units->year / thi->units->meter;
974: max[j] *= thi->units->year / thi->units->meter;
975: }
976: if (min[0] == 0.0) min[0] = 0.0;
977: PetscPrintf(comm,"|X|_2 %g %g <= u <= %g %g <= v <= %g %g <= c <= %g \n",(double)nrm2,(double)min[0],(double)max[0],(double)min[1],(double)max[1],(double)min[2],(double)max[2]);
978: {
979: PetscReal umin,umax,umean;
980: THISurfaceStatistics(dm,X,&umin,&umax,&umean);
981: umin *= thi->units->year / thi->units->meter;
982: umax *= thi->units->year / thi->units->meter;
983: umean *= thi->units->year / thi->units->meter;
984: PetscPrintf(comm,"Surface statistics: u in [%12.6e, %12.6e] mean %12.6e\n",(double)umin,(double)umax,(double)umean);
985: }
986: /* These values stay nondimensional */
987: PetscPrintf(comm,"Global eta range %g to %g converged range %g to %g\n",(double)thi->eta.min,(double)thi->eta.max,(double)thi->eta.cmin,(double)thi->eta.cmax);
988: PetscPrintf(comm,"Global beta2 range %g to %g converged range %g to %g\n",(double)thi->beta2.min,(double)thi->beta2.max,(double)thi->beta2.cmin,(double)thi->beta2.cmax);
989: }
990: return(0);
991: }
995: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info,Node **x,Mat J,Mat B,THI thi)
996: {
997: PetscInt xs,ys,xm,ym,i,j,q,l,ll;
998: PetscReal hx,hy;
999: PrmNode **prm;
1003: xs = info->ys;
1004: ys = info->xs;
1005: xm = info->ym;
1006: ym = info->xm;
1007: hx = thi->Lx / info->my;
1008: hy = thi->Ly / info->mx;
1010: MatZeroEntries(B);
1011: THIDAGetPrm(info->da,&prm);
1013: for (i=xs; i<xs+xm; i++) {
1014: for (j=ys; j<ys+ym; j++) {
1015: Node n[4];
1016: PrmNode pn[4];
1017: PetscScalar Ke[4*2][4*2];
1018: QuadExtract(prm,i,j,pn);
1019: QuadExtract(x,i,j,n);
1020: PetscMemzero(Ke,sizeof(Ke));
1021: for (q=0; q<4; q++) {
1022: PetscReal phi[4],dphi[4][2],jw,eta,deta,beta2,dbeta2;
1023: PetscScalar u,v,du[2],dv[2],h = 0,rbeta2 = 0;
1024: for (l=0; l<4; l++) {
1025: phi[l] = QuadQInterp[q][l];
1026: dphi[l][0] = QuadQDeriv[q][l][0]*2./hx;
1027: dphi[l][1] = QuadQDeriv[q][l][1]*2./hy;
1028: h += phi[l] * pn[l].h;
1029: rbeta2 += phi[l] * pn[l].beta2;
1030: }
1031: jw = 0.25*hx*hy / thi->rhog; /* rhog is only scaling */
1032: PointwiseNonlinearity2D(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
1033: THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
1034: for (l=0; l<4; l++) {
1035: const PetscReal pp = phi[l],*dp = dphi[l];
1036: for (ll=0; ll<4; ll++) {
1037: const PetscReal ppl = phi[ll],*dpl = dphi[ll];
1038: PetscScalar dgdu,dgdv;
1039: dgdu = 2.*du[0]*dpl[0] + dv[1]*dpl[0] + 0.5*(du[1]+dv[0])*dpl[1];
1040: dgdv = 2.*dv[1]*dpl[1] + du[0]*dpl[1] + 0.5*(du[1]+dv[0])*dpl[0];
1041: /* Picard part */
1042: Ke[l*2+0][ll*2+0] += dp[0]*jw*eta*4.*dpl[0] + dp[1]*jw*eta*dpl[1] + pp*jw*(beta2/h)*ppl*thi->ssa_friction_scale;
1043: Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0];
1044: Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1];
1045: Ke[l*2+1][ll*2+1] += dp[1]*jw*eta*4.*dpl[1] + dp[0]*jw*eta*dpl[0] + pp*jw*(beta2/h)*ppl*thi->ssa_friction_scale;
1046: /* extra Newton terms */
1047: 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]) + pp*jw*(dbeta2/h)*u*u*ppl*thi->ssa_friction_scale;
1048: 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]) + pp*jw*(dbeta2/h)*u*v*ppl*thi->ssa_friction_scale;
1049: 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]) + pp*jw*(dbeta2/h)*v*u*ppl*thi->ssa_friction_scale;
1050: 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]) + pp*jw*(dbeta2/h)*v*v*ppl*thi->ssa_friction_scale;
1051: }
1052: }
1053: }
1054: {
1055: const MatStencil rc[4] = {{0,i,j,0},{0,i+1,j,0},{0,i+1,j+1,0},{0,i,j+1,0}};
1056: MatSetValuesBlockedStencil(B,4,rc,4,rc,&Ke[0][0],ADD_VALUES);
1057: }
1058: }
1059: }
1060: THIDARestorePrm(info->da,&prm);
1062: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1063: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1064: MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1065: if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1066: return(0);
1067: }
1071: static PetscErrorCode THIJacobianLocal_3D(DMDALocalInfo *info,Node ***x,Mat B,THI thi,THIAssemblyMode amode)
1072: {
1073: PetscInt xs,ys,xm,ym,zm,i,j,k,q,l,ll;
1074: PetscReal hx,hy;
1075: PrmNode **prm;
1079: xs = info->zs;
1080: ys = info->ys;
1081: xm = info->zm;
1082: ym = info->ym;
1083: zm = info->xm;
1084: hx = thi->Lx / info->mz;
1085: hy = thi->Ly / info->my;
1087: MatZeroEntries(B);
1088: MatSetOption(B,MAT_SUBSET_OFF_PROC_ENTRIES,PETSC_TRUE);
1089: THIDAGetPrm(info->da,&prm);
1091: for (i=xs; i<xs+xm; i++) {
1092: for (j=ys; j<ys+ym; j++) {
1093: PrmNode pn[4];
1094: QuadExtract(prm,i,j,pn);
1095: for (k=0; k<zm-1; k++) {
1096: Node n[8];
1097: PetscReal zn[8],etabase = 0;
1098: PetscScalar Ke[8*2][8*2];
1099: PetscInt ls = 0;
1101: PrmHexGetZ(pn,k,zm,zn);
1102: HexExtract(x,i,j,k,n);
1103: PetscMemzero(Ke,sizeof(Ke));
1104: if (thi->no_slip && k == 0) {
1105: for (l=0; l<4; l++) n[l].u = n[l].v = 0;
1106: ls = 4;
1107: }
1108: for (q=0; q<8; q++) {
1109: PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta;
1110: PetscScalar du[3],dv[3],u,v;
1111: HexGrad(HexQDeriv[q],zn,dz);
1112: HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
1113: PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
1114: jw /= thi->rhog; /* residuals are scaled by this factor */
1115: if (q == 0) etabase = eta;
1116: for (l=ls; l<8; l++) { /* test functions */
1117: const PetscReal *PETSC_RESTRICT dp = dphi[l];
1118: #if USE_SSE2_KERNELS
1119: /* gcc (up to my 4.5 snapshot) is really bad at hoisting intrinsics so we do it manually */
1120: __m128d
1121: p4 = _mm_set1_pd(4),p2 = _mm_set1_pd(2),p05 = _mm_set1_pd(0.5),
1122: p42 = _mm_setr_pd(4,2),p24 = _mm_shuffle_pd(p42,p42,_MM_SHUFFLE2(0,1)),
1123: du0 = _mm_set1_pd(du[0]),du1 = _mm_set1_pd(du[1]),du2 = _mm_set1_pd(du[2]),
1124: dv0 = _mm_set1_pd(dv[0]),dv1 = _mm_set1_pd(dv[1]),dv2 = _mm_set1_pd(dv[2]),
1125: jweta = _mm_set1_pd(jw*eta),jwdeta = _mm_set1_pd(jw*deta),
1126: dp0 = _mm_set1_pd(dp[0]),dp1 = _mm_set1_pd(dp[1]),dp2 = _mm_set1_pd(dp[2]),
1127: dp0jweta = _mm_mul_pd(dp0,jweta),dp1jweta = _mm_mul_pd(dp1,jweta),dp2jweta = _mm_mul_pd(dp2,jweta),
1128: p4du0p2dv1 = _mm_add_pd(_mm_mul_pd(p4,du0),_mm_mul_pd(p2,dv1)), /* 4 du0 + 2 dv1 */
1129: p4dv1p2du0 = _mm_add_pd(_mm_mul_pd(p4,dv1),_mm_mul_pd(p2,du0)), /* 4 dv1 + 2 du0 */
1130: pdu2dv2 = _mm_unpacklo_pd(du2,dv2), /* [du2, dv2] */
1131: du1pdv0 = _mm_add_pd(du1,dv0), /* du1 + dv0 */
1132: t1 = _mm_mul_pd(dp0,p4du0p2dv1), /* dp0 (4 du0 + 2 dv1) */
1133: t2 = _mm_mul_pd(dp1,p4dv1p2du0); /* dp1 (4 dv1 + 2 du0) */
1135: #endif
1136: #if defined COMPUTE_LOWER_TRIANGULAR /* The element matrices are always symmetric so computing the lower-triangular part is not necessary */
1137: for (ll=ls; ll<8; ll++) { /* trial functions */
1138: #else
1139: for (ll=l; ll<8; ll++) {
1140: #endif
1141: const PetscReal *PETSC_RESTRICT dpl = dphi[ll];
1142: if (amode == THIASSEMBLY_TRIDIAGONAL && (l-ll)%4) continue; /* these entries would not be inserted */
1143: #if !USE_SSE2_KERNELS
1144: /* The analytic Jacobian in nice, easy-to-read form */
1145: {
1146: PetscScalar dgdu,dgdv;
1147: dgdu = 2.*du[0]*dpl[0] + dv[1]*dpl[0] + 0.5*(du[1]+dv[0])*dpl[1] + 0.5*du[2]*dpl[2];
1148: dgdv = 2.*dv[1]*dpl[1] + du[0]*dpl[1] + 0.5*(du[1]+dv[0])*dpl[0] + 0.5*dv[2]*dpl[2];
1149: /* Picard part */
1150: 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];
1151: Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0];
1152: Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1];
1153: 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];
1154: /* extra Newton terms */
1155: 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];
1156: 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];
1157: 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];
1158: 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];
1159: }
1160: #else
1161: /* This SSE2 code is an exact replica of above, but uses explicit packed instructions for some speed
1162: * benefit. On my hardware, these intrinsics are almost twice as fast as above, reducing total assembly cost
1163: * by 25 to 30 percent. */
1164: {
1165: __m128d
1166: keu = _mm_loadu_pd(&Ke[l*2+0][ll*2+0]),
1167: kev = _mm_loadu_pd(&Ke[l*2+1][ll*2+0]),
1168: dpl01 = _mm_loadu_pd(&dpl[0]),dpl10 = _mm_shuffle_pd(dpl01,dpl01,_MM_SHUFFLE2(0,1)),dpl2 = _mm_set_sd(dpl[2]),
1169: t0,t3,pdgduv;
1170: keu = _mm_add_pd(keu,_mm_add_pd(_mm_mul_pd(_mm_mul_pd(dp0jweta,p42),dpl01),
1171: _mm_add_pd(_mm_mul_pd(dp1jweta,dpl10),
1172: _mm_mul_pd(dp2jweta,dpl2))));
1173: kev = _mm_add_pd(kev,_mm_add_pd(_mm_mul_pd(_mm_mul_pd(dp1jweta,p24),dpl01),
1174: _mm_add_pd(_mm_mul_pd(dp0jweta,dpl10),
1175: _mm_mul_pd(dp2jweta,_mm_shuffle_pd(dpl2,dpl2,_MM_SHUFFLE2(0,1))))));
1176: pdgduv = _mm_mul_pd(p05,_mm_add_pd(_mm_add_pd(_mm_mul_pd(p42,_mm_mul_pd(du0,dpl01)),
1177: _mm_mul_pd(p24,_mm_mul_pd(dv1,dpl01))),
1178: _mm_add_pd(_mm_mul_pd(du1pdv0,dpl10),
1179: _mm_mul_pd(pdu2dv2,_mm_set1_pd(dpl[2]))))); /* [dgdu, dgdv] */
1180: t0 = _mm_mul_pd(jwdeta,pdgduv); /* jw deta [dgdu, dgdv] */
1181: t3 = _mm_mul_pd(t0,du1pdv0); /* t0 (du1 + dv0) */
1182: _mm_storeu_pd(&Ke[l*2+0][ll*2+0],_mm_add_pd(keu,_mm_add_pd(_mm_mul_pd(t1,t0),
1183: _mm_add_pd(_mm_mul_pd(dp1,t3),
1184: _mm_mul_pd(t0,_mm_mul_pd(dp2,du2))))));
1185: _mm_storeu_pd(&Ke[l*2+1][ll*2+0],_mm_add_pd(kev,_mm_add_pd(_mm_mul_pd(t2,t0),
1186: _mm_add_pd(_mm_mul_pd(dp0,t3),
1187: _mm_mul_pd(t0,_mm_mul_pd(dp2,dv2))))));
1188: }
1189: #endif
1190: }
1191: }
1192: }
1193: if (k == 0) { /* on a bottom face */
1194: if (thi->no_slip) {
1195: const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1);
1196: 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);
1197: Ke[0][0] = thi->dirichlet_scale*diagu;
1198: Ke[1][1] = thi->dirichlet_scale*diagv;
1199: } else {
1200: for (q=0; q<4; q++) {
1201: const PetscReal jw = 0.25*hx*hy/thi->rhog,*phi = QuadQInterp[q];
1202: PetscScalar u =0,v=0,rbeta2=0;
1203: PetscReal beta2,dbeta2;
1204: for (l=0; l<4; l++) {
1205: u += phi[l]*n[l].u;
1206: v += phi[l]*n[l].v;
1207: rbeta2 += phi[l]*pn[l].beta2;
1208: }
1209: THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
1210: for (l=0; l<4; l++) {
1211: const PetscReal pp = phi[l];
1212: for (ll=0; ll<4; ll++) {
1213: const PetscReal ppl = phi[ll];
1214: Ke[l*2+0][ll*2+0] += pp*jw*beta2*ppl + pp*jw*dbeta2*u*u*ppl;
1215: Ke[l*2+0][ll*2+1] += pp*jw*dbeta2*u*v*ppl;
1216: Ke[l*2+1][ll*2+0] += pp*jw*dbeta2*v*u*ppl;
1217: Ke[l*2+1][ll*2+1] += pp*jw*beta2*ppl + pp*jw*dbeta2*v*v*ppl;
1218: }
1219: }
1220: }
1221: }
1222: }
1223: {
1224: const MatStencil rc[8] = {{i,j,k,0},{i+1,j,k,0},{i+1,j+1,k,0},{i,j+1,k,0},{i,j,k+1,0},{i+1,j,k+1,0},{i+1,j+1,k+1,0},{i,j+1,k+1,0}};
1225: if (amode == THIASSEMBLY_TRIDIAGONAL) {
1226: for (l=0; l<4; l++) { /* Copy out each of the blocks, discarding horizontal coupling */
1227: const PetscInt l4 = l+4;
1228: const MatStencil rcl[2] = {{rc[l].k,rc[l].j,rc[l].i,0},{rc[l4].k,rc[l4].j,rc[l4].i,0}};
1229: #if defined COMPUTE_LOWER_TRIANGULAR
1230: const PetscScalar Kel[4][4] = {{Ke[2*l+0][2*l+0] ,Ke[2*l+0][2*l+1] ,Ke[2*l+0][2*l4+0] ,Ke[2*l+0][2*l4+1]},
1231: {Ke[2*l+1][2*l+0] ,Ke[2*l+1][2*l+1] ,Ke[2*l+1][2*l4+0] ,Ke[2*l+1][2*l4+1]},
1232: {Ke[2*l4+0][2*l+0],Ke[2*l4+0][2*l+1],Ke[2*l4+0][2*l4+0],Ke[2*l4+0][2*l4+1]},
1233: {Ke[2*l4+1][2*l+0],Ke[2*l4+1][2*l+1],Ke[2*l4+1][2*l4+0],Ke[2*l4+1][2*l4+1]}};
1234: #else
1235: /* Same as above except for the lower-left block */
1236: const PetscScalar Kel[4][4] = {{Ke[2*l+0][2*l+0] ,Ke[2*l+0][2*l+1] ,Ke[2*l+0][2*l4+0] ,Ke[2*l+0][2*l4+1]},
1237: {Ke[2*l+1][2*l+0] ,Ke[2*l+1][2*l+1] ,Ke[2*l+1][2*l4+0] ,Ke[2*l+1][2*l4+1]},
1238: {Ke[2*l+0][2*l4+0],Ke[2*l+1][2*l4+0],Ke[2*l4+0][2*l4+0],Ke[2*l4+0][2*l4+1]},
1239: {Ke[2*l+0][2*l4+1],Ke[2*l+1][2*l4+1],Ke[2*l4+1][2*l4+0],Ke[2*l4+1][2*l4+1]}};
1240: #endif
1241: MatSetValuesBlockedStencil(B,2,rcl,2,rcl,&Kel[0][0],ADD_VALUES);
1242: }
1243: } else {
1244: #if !defined COMPUTE_LOWER_TRIANGULAR /* fill in lower-triangular part, this is really cheap compared to computing the entries */
1245: for (l=0; l<8; l++) {
1246: for (ll=l+1; ll<8; ll++) {
1247: Ke[ll*2+0][l*2+0] = Ke[l*2+0][ll*2+0];
1248: Ke[ll*2+1][l*2+0] = Ke[l*2+0][ll*2+1];
1249: Ke[ll*2+0][l*2+1] = Ke[l*2+1][ll*2+0];
1250: Ke[ll*2+1][l*2+1] = Ke[l*2+1][ll*2+1];
1251: }
1252: }
1253: #endif
1254: MatSetValuesBlockedStencil(B,8,rc,8,rc,&Ke[0][0],ADD_VALUES);
1255: }
1256: }
1257: }
1258: }
1259: }
1260: THIDARestorePrm(info->da,&prm);
1262: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1263: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1264: MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1265: if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1266: return(0);
1267: }
1271: static PetscErrorCode THIJacobianLocal_3D_Full(DMDALocalInfo *info,Node ***x,Mat A,Mat B,THI thi)
1272: {
1276: THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_FULL);
1277: return(0);
1278: }
1282: static PetscErrorCode THIJacobianLocal_3D_Tridiagonal(DMDALocalInfo *info,Node ***x,Mat A,Mat B,THI thi)
1283: {
1287: THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_TRIDIAGONAL);
1288: return(0);
1289: }
1293: static PetscErrorCode DMRefineHierarchy_THI(DM dac0,PetscInt nlevels,DM hierarchy[])
1294: {
1295: PetscErrorCode ierr;
1296: THI thi;
1297: PetscInt dim,M,N,m,n,s,dof;
1298: DM dac,daf;
1299: DMDAStencilType st;
1300: DM_DA *ddf,*ddc;
1303: PetscObjectQuery((PetscObject)dac0,"THI",(PetscObject*)&thi);
1304: if (!thi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot refine this DMDA, missing composed THI instance");
1305: if (nlevels > 1) {
1306: DMRefineHierarchy(dac0,nlevels-1,hierarchy);
1307: dac = hierarchy[nlevels-2];
1308: } else {
1309: dac = dac0;
1310: }
1311: DMDAGetInfo(dac,&dim, &N,&M,0, &n,&m,0, &dof,&s,0,0,0,&st);
1312: if (dim != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"This function can only refine 2D DMDAs");
1314: /* Creates a 3D DMDA with the same map-plane layout as the 2D one, with contiguous columns */
1315: DMDACreate3d(PetscObjectComm((PetscObject)dac),DM_BOUNDARY_NONE,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,thi->zlevels,N,M,1,n,m,dof,s,NULL,NULL,NULL,&daf);
1317: daf->ops->creatematrix = dac->ops->creatematrix;
1318: daf->ops->createinterpolation = dac->ops->createinterpolation;
1319: daf->ops->getcoloring = dac->ops->getcoloring;
1320: ddf = (DM_DA*)daf->data;
1321: ddc = (DM_DA*)dac->data;
1322: ddf->interptype = ddc->interptype;
1324: DMDASetFieldName(daf,0,"x-velocity");
1325: DMDASetFieldName(daf,1,"y-velocity");
1327: hierarchy[nlevels-1] = daf;
1328: return(0);
1329: }
1333: static PetscErrorCode DMCreateInterpolation_DA_THI(DM dac,DM daf,Mat *A,Vec *scale)
1334: {
1336: PetscInt dim;
1343: DMDAGetInfo(daf,&dim,0,0,0,0,0,0,0,0,0,0,0,0);
1344: if (dim == 2) {
1345: /* We are in the 2D problem and use normal DMDA interpolation */
1346: DMCreateInterpolation(dac,daf,A,scale);
1347: } else {
1348: PetscInt i,j,k,xs,ys,zs,xm,ym,zm,mx,my,mz,rstart,cstart;
1349: Mat B;
1351: DMDAGetInfo(daf,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
1352: DMDAGetCorners(daf,&zs,&ys,&xs,&zm,&ym,&xm);
1353: if (zs != 0) SETERRQ(PETSC_COMM_SELF,1,"unexpected");
1354: MatCreate(PetscObjectComm((PetscObject)daf),&B);
1355: MatSetSizes(B,xm*ym*zm,xm*ym,mx*my*mz,mx*my);
1357: MatSetType(B,MATAIJ);
1358: MatSeqAIJSetPreallocation(B,1,NULL);
1359: MatMPIAIJSetPreallocation(B,1,NULL,0,NULL);
1360: MatGetOwnershipRange(B,&rstart,NULL);
1361: MatGetOwnershipRangeColumn(B,&cstart,NULL);
1362: for (i=xs; i<xs+xm; i++) {
1363: for (j=ys; j<ys+ym; j++) {
1364: for (k=zs; k<zs+zm; k++) {
1365: PetscInt i2 = i*ym+j,i3 = i2*zm+k;
1366: PetscScalar val = ((k == 0 || k == mz-1) ? 0.5 : 1.) / (mz-1.); /* Integration using trapezoid rule */
1367: MatSetValue(B,cstart+i3,rstart+i2,val,INSERT_VALUES);
1368: }
1369: }
1370: }
1371: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1372: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1373: MatCreateMAIJ(B,sizeof(Node)/sizeof(PetscScalar),A);
1374: MatDestroy(&B);
1375: }
1376: return(0);
1377: }
1381: static PetscErrorCode DMCreateMatrix_THI_Tridiagonal(DM da,Mat *J)
1382: {
1383: PetscErrorCode ierr;
1384: Mat A;
1385: PetscInt xm,ym,zm,dim,dof = 2,starts[3],dims[3];
1386: ISLocalToGlobalMapping ltog;
1389: DMDAGetInfo(da,&dim, 0,0,0, 0,0,0, 0,0,0,0,0,0);
1390: if (dim != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected DMDA to be 3D");
1391: DMDAGetCorners(da,0,0,0,&zm,&ym,&xm);
1392: DMGetLocalToGlobalMapping(da,<og);
1393: MatCreate(PetscObjectComm((PetscObject)da),&A);
1394: MatSetSizes(A,dof*xm*ym*zm,dof*xm*ym*zm,PETSC_DETERMINE,PETSC_DETERMINE);
1395: MatSetType(A,da->mattype);
1396: MatSetFromOptions(A);
1397: MatSeqAIJSetPreallocation(A,3*2,NULL);
1398: MatMPIAIJSetPreallocation(A,3*2,NULL,0,NULL);
1399: MatSeqBAIJSetPreallocation(A,2,3,NULL);
1400: MatMPIBAIJSetPreallocation(A,2,3,NULL,0,NULL);
1401: MatSeqSBAIJSetPreallocation(A,2,2,NULL);
1402: MatMPISBAIJSetPreallocation(A,2,2,NULL,0,NULL);
1403: MatSetLocalToGlobalMapping(A,ltog,ltog);
1404: DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);
1405: MatSetStencil(A,dim,dims,starts,dof);
1406: *J = A;
1407: return(0);
1408: }
1412: static PetscErrorCode THIDAVecView_VTK_XML(THI thi,DM da,Vec X,const char filename[])
1413: {
1414: const PetscInt dof = 2;
1415: Units units = thi->units;
1416: MPI_Comm comm;
1417: PetscErrorCode ierr;
1418: PetscViewer viewer;
1419: PetscMPIInt rank,size,tag,nn,nmax;
1420: PetscInt mx,my,mz,r,range[6];
1421: const PetscScalar *x;
1424: PetscObjectGetComm((PetscObject)thi,&comm);
1425: DMDAGetInfo(da,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
1426: MPI_Comm_size(comm,&size);
1427: MPI_Comm_rank(comm,&rank);
1428: PetscViewerASCIIOpen(comm,filename,&viewer);
1429: PetscViewerASCIIPrintf(viewer,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");
1430: PetscViewerASCIIPrintf(viewer," <StructuredGrid WholeExtent=\"%d %D %d %D %d %D\">\n",0,mz-1,0,my-1,0,mx-1);
1432: DMDAGetCorners(da,range,range+1,range+2,range+3,range+4,range+5);
1433: PetscMPIIntCast(range[3]*range[4]*range[5]*dof,&nn);
1434: MPI_Reduce(&nn,&nmax,1,MPI_INT,MPI_MAX,0,comm);
1435: tag = ((PetscObject) viewer)->tag;
1436: VecGetArrayRead(X,&x);
1437: if (!rank) {
1438: PetscScalar *array;
1439: PetscMalloc1(nmax,&array);
1440: for (r=0; r<size; r++) {
1441: PetscInt i,j,k,xs,xm,ys,ym,zs,zm;
1442: const PetscScalar *ptr;
1443: MPI_Status status;
1444: if (r) {
1445: MPI_Recv(range,6,MPIU_INT,r,tag,comm,MPI_STATUS_IGNORE);
1446: }
1447: zs = range[0];ys = range[1];xs = range[2];zm = range[3];ym = range[4];xm = range[5];
1448: if (xm*ym*zm*dof > nmax) SETERRQ(PETSC_COMM_SELF,1,"should not happen");
1449: if (r) {
1450: MPI_Recv(array,nmax,MPIU_SCALAR,r,tag,comm,&status);
1451: MPI_Get_count(&status,MPIU_SCALAR,&nn);
1452: if (nn != xm*ym*zm*dof) SETERRQ(PETSC_COMM_SELF,1,"should not happen");
1453: ptr = array;
1454: } else ptr = x;
1455: PetscViewerASCIIPrintf(viewer," <Piece Extent=\"%D %D %D %D %D %D\">\n",zs,zs+zm-1,ys,ys+ym-1,xs,xs+xm-1);
1457: PetscViewerASCIIPrintf(viewer," <Points>\n");
1458: PetscViewerASCIIPrintf(viewer," <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n");
1459: for (i=xs; i<xs+xm; i++) {
1460: for (j=ys; j<ys+ym; j++) {
1461: for (k=zs; k<zs+zm; k++) {
1462: PrmNode p;
1463: PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my,zz;
1464: thi->initialize(thi,xx,yy,&p);
1465: zz = PetscRealPart(p.b) + PetscRealPart(p.h)*k/(mz-1);
1466: PetscViewerASCIIPrintf(viewer,"%f %f %f\n",(double)xx,(double)yy,(double)zz);
1467: }
1468: }
1469: }
1470: PetscViewerASCIIPrintf(viewer," </DataArray>\n");
1471: PetscViewerASCIIPrintf(viewer," </Points>\n");
1473: PetscViewerASCIIPrintf(viewer," <PointData>\n");
1474: PetscViewerASCIIPrintf(viewer," <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n");
1475: for (i=0; i<nn; i+=dof) {
1476: PetscViewerASCIIPrintf(viewer,"%f %f %f\n",(double)(PetscRealPart(ptr[i])*units->year/units->meter),(double)(PetscRealPart(ptr[i+1])*units->year/units->meter),0.0);
1477: }
1478: PetscViewerASCIIPrintf(viewer," </DataArray>\n");
1480: PetscViewerASCIIPrintf(viewer," <DataArray type=\"Int32\" Name=\"rank\" NumberOfComponents=\"1\" format=\"ascii\">\n");
1481: for (i=0; i<nn; i+=dof) {
1482: PetscViewerASCIIPrintf(viewer,"%D\n",r);
1483: }
1484: PetscViewerASCIIPrintf(viewer," </DataArray>\n");
1485: PetscViewerASCIIPrintf(viewer," </PointData>\n");
1487: PetscViewerASCIIPrintf(viewer," </Piece>\n");
1488: }
1489: PetscFree(array);
1490: } else {
1491: MPI_Send(range,6,MPIU_INT,0,tag,comm);
1492: MPI_Send((PetscScalar*)x,nn,MPIU_SCALAR,0,tag,comm);
1493: }
1494: VecRestoreArrayRead(X,&x);
1495: PetscViewerASCIIPrintf(viewer," </StructuredGrid>\n");
1496: PetscViewerASCIIPrintf(viewer,"</VTKFile>\n");
1497: PetscViewerDestroy(&viewer);
1498: return(0);
1499: }
1503: int main(int argc,char *argv[])
1504: {
1505: MPI_Comm comm;
1506: THI thi;
1508: DM da;
1509: SNES snes;
1511: PetscInitialize(&argc,&argv,0,help);
1512: comm = PETSC_COMM_WORLD;
1514: THICreate(comm,&thi);
1515: {
1516: PetscInt M = 3,N = 3,P = 2;
1517: PetscOptionsBegin(comm,NULL,"Grid resolution options","");
1518: {
1519: PetscOptionsInt("-M","Number of elements in x-direction on coarse level","",M,&M,NULL);
1520: N = M;
1521: PetscOptionsInt("-N","Number of elements in y-direction on coarse level (if different from M)","",N,&N,NULL);
1522: if (thi->coarse2d) {
1523: PetscOptionsInt("-zlevels","Number of elements in z-direction on fine level","",thi->zlevels,&thi->zlevels,NULL);
1524: } else {
1525: PetscOptionsInt("-P","Number of elements in z-direction on coarse level","",P,&P,NULL);
1526: }
1527: }
1528: PetscOptionsEnd();
1529: if (thi->coarse2d) {
1530: DMDACreate2d(comm,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX,-N,-M,PETSC_DETERMINE,PETSC_DETERMINE,sizeof(Node)/sizeof(PetscScalar),1,0,0,&da);
1532: da->ops->refinehierarchy = DMRefineHierarchy_THI;
1533: da->ops->createinterpolation = DMCreateInterpolation_DA_THI;
1535: PetscObjectCompose((PetscObject)da,"THI",(PetscObject)thi);
1536: } else {
1537: 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);
1538: }
1539: DMDASetFieldName(da,0,"x-velocity");
1540: DMDASetFieldName(da,1,"y-velocity");
1541: }
1542: THISetUpDM(thi,da);
1543: if (thi->tridiagonal) da->ops->creatematrix = DMCreateMatrix_THI_Tridiagonal;
1545: { /* Set the fine level matrix type if -da_refine */
1546: PetscInt rlevel,clevel;
1547: DMGetRefineLevel(da,&rlevel);
1548: DMGetCoarsenLevel(da,&clevel);
1549: if (rlevel - clevel > 0) {DMSetMatType(da,thi->mattype);}
1550: }
1552: DMDASNESSetFunctionLocal(da,ADD_VALUES,(DMDASNESFunction)THIFunctionLocal,thi);
1553: if (thi->tridiagonal) {
1554: DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)THIJacobianLocal_3D_Tridiagonal,thi);
1555: } else {
1556: DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)THIJacobianLocal_3D_Full,thi);
1557: }
1558: DMCoarsenHookAdd(da,DMCoarsenHook_THI,NULL,thi);
1559: DMRefineHookAdd(da,DMRefineHook_THI,NULL,thi);
1561: DMSetApplicationContext(da,thi);
1563: SNESCreate(comm,&snes);
1564: SNESSetDM(snes,da);
1565: DMDestroy(&da);
1566: SNESSetComputeInitialGuess(snes,THIInitial,NULL);
1567: SNESSetFromOptions(snes);
1569: SNESSolve(snes,NULL,NULL);
1571: THISolveStatistics(thi,snes,0,"Full");
1573: {
1574: PetscBool flg;
1575: char filename[PETSC_MAX_PATH_LEN] = "";
1576: PetscOptionsGetString(NULL,NULL,"-o",filename,sizeof(filename),&flg);
1577: if (flg) {
1578: Vec X;
1579: DM dm;
1580: SNESGetSolution(snes,&X);
1581: SNESGetDM(snes,&dm);
1582: THIDAVecView_VTK_XML(thi,dm,X,filename);
1583: }
1584: }
1586: DMDestroy(&da);
1587: SNESDestroy(&snes);
1588: THIDestroy(&thi);
1589: PetscFinalize();
1590: return 0;
1591: }