Actual source code: ex48.c
petsc-3.8.4 2018-03-24
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) || defined (PETSC_HAVE_WINDOWS_COMPILERS) || defined (__PGI)
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: #if !defined NO_SSE2 \
79: && !defined PETSC_USE_COMPLEX \
80: && !defined PETSC_USE_REAL_SINGLE \
81: && !defined PETSC_USE_REAL___FLOAT128 \
82: && !defined PETSC_USE_REAL___FP16 \
83: && defined __SSE2__
84: #define USE_SSE2_KERNELS 1
85: #else
86: #define USE_SSE2_KERNELS 0
87: #endif
89: static PetscClassId THI_CLASSID;
91: typedef enum {QUAD_GAUSS,QUAD_LOBATTO} QuadratureType;
92: static const char *QuadratureTypes[] = {"gauss","lobatto","QuadratureType","QUAD_",0};
93: PETSC_UNUSED static const PetscReal HexQWeights[8] = {1,1,1,1,1,1,1,1};
94: PETSC_UNUSED static const PetscReal HexQNodes[] = {-0.57735026918962573, 0.57735026918962573};
95: #define G 0.57735026918962573
96: #define H (0.5*(1.+G))
97: #define L (0.5*(1.-G))
98: #define M (-0.5)
99: #define P (0.5)
100: /* Special quadrature: Lobatto in horizontal, Gauss in vertical */
101: static const PetscReal HexQInterp_Lobatto[8][8] = {{H,0,0,0,L,0,0,0},
102: {0,H,0,0,0,L,0,0},
103: {0,0,H,0,0,0,L,0},
104: {0,0,0,H,0,0,0,L},
105: {L,0,0,0,H,0,0,0},
106: {0,L,0,0,0,H,0,0},
107: {0,0,L,0,0,0,H,0},
108: {0,0,0,L,0,0,0,H}};
109: static const PetscReal HexQDeriv_Lobatto[8][8][3] = {
110: {{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} },
111: {{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} },
112: {{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} },
113: {{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}},
114: {{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} },
115: {{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} },
116: {{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} },
117: {{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}}};
118: /* Stanndard Gauss */
119: 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},
120: {L*H*H,H*H*H,H*L*H,L*L*H, L*H*L,H*H*L,H*L*L,L*L*L},
121: {L*L*H,H*L*H,H*H*H,L*H*H, L*L*L,H*L*L,H*H*L,L*H*L},
122: {H*L*H,L*L*H,L*H*H,H*H*H, H*L*L,L*L*L,L*H*L,H*H*L},
123: {H*H*L,L*H*L,L*L*L,H*L*L, H*H*H,L*H*H,L*L*H,H*L*H},
124: {L*H*L,H*H*L,H*L*L,L*L*L, L*H*H,H*H*H,H*L*H,L*L*H},
125: {L*L*L,H*L*L,H*H*L,L*H*L, L*L*H,H*L*H,H*H*H,L*H*H},
126: {H*L*L,L*L*L,L*H*L,H*H*L, H*L*H,L*L*H,L*H*H,H*H*H}};
127: static const PetscReal HexQDeriv_Gauss[8][8][3] = {
128: {{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}},
129: {{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}},
130: {{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}},
131: {{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}},
132: {{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}},
133: {{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}},
134: {{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}},
135: {{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}}};
136: static const PetscReal (*HexQInterp)[8],(*HexQDeriv)[8][3];
137: /* Standard 2x2 Gauss quadrature for the bottom layer. */
138: static const PetscReal QuadQInterp[4][4] = {{H*H,L*H,L*L,H*L},
139: {L*H,H*H,H*L,L*L},
140: {L*L,H*L,H*H,L*H},
141: {H*L,L*L,L*H,H*H}};
142: static const PetscReal QuadQDeriv[4][4][2] = {
143: {{M*H,M*H},{P*H,M*L},{P*L,P*L},{M*L,P*H}},
144: {{M*H,M*L},{P*H,M*H},{P*L,P*H},{M*L,P*L}},
145: {{M*L,M*L},{P*L,M*H},{P*H,P*H},{M*H,P*L}},
146: {{M*L,M*H},{P*L,M*L},{P*H,P*L},{M*H,P*H}}};
147: #undef G
148: #undef H
149: #undef L
150: #undef M
151: #undef P
153: #define HexExtract(x,i,j,k,n) do { \
154: (n)[0] = (x)[i][j][k]; \
155: (n)[1] = (x)[i+1][j][k]; \
156: (n)[2] = (x)[i+1][j+1][k]; \
157: (n)[3] = (x)[i][j+1][k]; \
158: (n)[4] = (x)[i][j][k+1]; \
159: (n)[5] = (x)[i+1][j][k+1]; \
160: (n)[6] = (x)[i+1][j+1][k+1]; \
161: (n)[7] = (x)[i][j+1][k+1]; \
162: } while (0)
164: #define HexExtractRef(x,i,j,k,n) do { \
165: (n)[0] = &(x)[i][j][k]; \
166: (n)[1] = &(x)[i+1][j][k]; \
167: (n)[2] = &(x)[i+1][j+1][k]; \
168: (n)[3] = &(x)[i][j+1][k]; \
169: (n)[4] = &(x)[i][j][k+1]; \
170: (n)[5] = &(x)[i+1][j][k+1]; \
171: (n)[6] = &(x)[i+1][j+1][k+1]; \
172: (n)[7] = &(x)[i][j+1][k+1]; \
173: } while (0)
175: #define QuadExtract(x,i,j,n) do { \
176: (n)[0] = (x)[i][j]; \
177: (n)[1] = (x)[i+1][j]; \
178: (n)[2] = (x)[i+1][j+1]; \
179: (n)[3] = (x)[i][j+1]; \
180: } while (0)
182: static void HexGrad(const PetscReal dphi[][3],const PetscReal zn[],PetscReal dz[])
183: {
184: PetscInt i;
185: dz[0] = dz[1] = dz[2] = 0;
186: for (i=0; i<8; i++) {
187: dz[0] += dphi[i][0] * zn[i];
188: dz[1] += dphi[i][1] * zn[i];
189: dz[2] += dphi[i][2] * zn[i];
190: }
191: }
193: 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)
194: {
195: const PetscReal jac[3][3] = {{hx/2,0,0}, {0,hy/2,0}, {dz[0],dz[1],dz[2]}};
196: 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]}};
197: const PetscReal jdet = jac[0][0]*jac[1][1]*jac[2][2];
198: PetscInt i;
200: for (i=0; i<8; i++) {
201: const PetscReal *dphir = HexQDeriv[q][i];
202: phi[i] = HexQInterp[q][i];
203: dphi[i][0] = dphir[0]*ijac[0][0] + dphir[1]*ijac[1][0] + dphir[2]*ijac[2][0];
204: dphi[i][1] = dphir[0]*ijac[0][1] + dphir[1]*ijac[1][1] + dphir[2]*ijac[2][1];
205: dphi[i][2] = dphir[0]*ijac[0][2] + dphir[1]*ijac[1][2] + dphir[2]*ijac[2][2];
206: }
207: *jw = 1.0 * jdet;
208: }
210: typedef struct _p_THI *THI;
211: typedef struct _n_Units *Units;
213: typedef struct {
214: PetscScalar u,v;
215: } Node;
217: typedef struct {
218: PetscScalar b; /* bed */
219: PetscScalar h; /* thickness */
220: PetscScalar beta2; /* friction */
221: } PrmNode;
223: typedef struct {
224: PetscReal min,max,cmin,cmax;
225: } PRange;
227: typedef enum {THIASSEMBLY_TRIDIAGONAL,THIASSEMBLY_FULL} THIAssemblyMode;
229: struct _p_THI {
230: PETSCHEADER(int);
231: void (*initialize)(THI,PetscReal x,PetscReal y,PrmNode *p);
232: PetscInt zlevels;
233: PetscReal Lx,Ly,Lz; /* Model domain */
234: PetscReal alpha; /* Bed angle */
235: Units units;
236: PetscReal dirichlet_scale;
237: PetscReal ssa_friction_scale;
238: PRange eta;
239: PRange beta2;
240: struct {
241: PetscReal Bd2,eps,exponent;
242: } viscosity;
243: struct {
244: PetscReal irefgam,eps2,exponent,refvel,epsvel;
245: } friction;
246: PetscReal rhog;
247: PetscBool no_slip;
248: PetscBool tridiagonal;
249: PetscBool coarse2d;
250: PetscBool verbose;
251: MatType mattype;
252: };
254: struct _n_Units {
255: /* fundamental */
256: PetscReal meter;
257: PetscReal kilogram;
258: PetscReal second;
259: /* derived */
260: PetscReal Pascal;
261: PetscReal year;
262: };
264: static PetscErrorCode THIJacobianLocal_3D_Full(DMDALocalInfo*,Node***,Mat,Mat,THI);
265: static PetscErrorCode THIJacobianLocal_3D_Tridiagonal(DMDALocalInfo*,Node***,Mat,Mat,THI);
266: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo*,Node**,Mat,Mat,THI);
268: static void PrmHexGetZ(const PrmNode pn[],PetscInt k,PetscInt zm,PetscReal zn[])
269: {
270: const PetscScalar zm1 = zm-1,
271: znl[8] = {pn[0].b + pn[0].h*(PetscScalar)k/zm1,
272: pn[1].b + pn[1].h*(PetscScalar)k/zm1,
273: pn[2].b + pn[2].h*(PetscScalar)k/zm1,
274: pn[3].b + pn[3].h*(PetscScalar)k/zm1,
275: pn[0].b + pn[0].h*(PetscScalar)(k+1)/zm1,
276: pn[1].b + pn[1].h*(PetscScalar)(k+1)/zm1,
277: pn[2].b + pn[2].h*(PetscScalar)(k+1)/zm1,
278: pn[3].b + pn[3].h*(PetscScalar)(k+1)/zm1};
279: PetscInt i;
280: for (i=0; i<8; i++) zn[i] = PetscRealPart(znl[i]);
281: }
283: /* Tests A and C are from the ISMIP-HOM paper (Pattyn et al. 2008) */
284: static void THIInitialize_HOM_A(THI thi,PetscReal x,PetscReal y,PrmNode *p)
285: {
286: Units units = thi->units;
287: PetscReal s = -x*PetscSinReal(thi->alpha);
289: p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx) * PetscSinReal(y*2*PETSC_PI/thi->Ly);
290: p->h = s - p->b;
291: p->beta2 = 1e30;
292: }
294: static void THIInitialize_HOM_C(THI thi,PetscReal x,PetscReal y,PrmNode *p)
295: {
296: Units units = thi->units;
297: PetscReal s = -x*PetscSinReal(thi->alpha);
299: p->b = s - 1000*units->meter;
300: p->h = s - p->b;
301: /* tau_b = beta2 v is a stress (Pa) */
302: p->beta2 = 1000 * (1 + PetscSinReal(x*2*PETSC_PI/thi->Lx)*PetscSinReal(y*2*PETSC_PI/thi->Ly)) * units->Pascal * units->year / units->meter;
303: }
305: /* These are just toys */
307: /* Same bed as test A, free slip everywhere except for a discontinuous jump to a circular sticky region in the middle. */
308: static void THIInitialize_HOM_X(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
309: {
310: Units units = thi->units;
311: PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
312: PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
313: p->b = s - 1000*units->meter + 500*units->meter*PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
314: p->h = s - p->b;
315: p->beta2 = 1000 * (r < 1 ? 2 : 0) * units->Pascal * units->year / units->meter;
316: }
318: /* Like Z, but with 200 meter cliffs */
319: static void THIInitialize_HOM_Y(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
320: {
321: Units units = thi->units;
322: PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
323: PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
325: p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
326: if (PetscRealPart(p->b) > -700*units->meter) p->b += 200*units->meter;
327: p->h = s - p->b;
328: 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;
329: }
331: /* Same bed as A, smoothly varying slipperiness, similar to MATLAB's "sombrero" (uncorrelated with bathymetry) */
332: static void THIInitialize_HOM_Z(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
333: {
334: Units units = thi->units;
335: PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
336: PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
338: p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
339: p->h = s - p->b;
340: 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;
341: }
343: static void THIFriction(THI thi,PetscReal rbeta2,PetscReal gam,PetscReal *beta2,PetscReal *dbeta2)
344: {
345: if (thi->friction.irefgam == 0) {
346: Units units = thi->units;
347: thi->friction.irefgam = 1./(0.5*PetscSqr(thi->friction.refvel * units->meter / units->year));
348: thi->friction.eps2 = 0.5*PetscSqr(thi->friction.epsvel * units->meter / units->year) * thi->friction.irefgam;
349: }
350: if (thi->friction.exponent == 0) {
351: *beta2 = rbeta2;
352: *dbeta2 = 0;
353: } else {
354: *beta2 = rbeta2 * PetscPowReal(thi->friction.eps2 + gam*thi->friction.irefgam,thi->friction.exponent);
355: *dbeta2 = thi->friction.exponent * *beta2 / (thi->friction.eps2 + gam*thi->friction.irefgam) * thi->friction.irefgam;
356: }
357: }
359: static void THIViscosity(THI thi,PetscReal gam,PetscReal *eta,PetscReal *deta)
360: {
361: PetscReal Bd2,eps,exponent;
362: if (thi->viscosity.Bd2 == 0) {
363: Units units = thi->units;
364: const PetscReal
365: n = 3., /* Glen exponent */
366: p = 1. + 1./n, /* for Stokes */
367: A = 1.e-16 * PetscPowReal(units->Pascal,-n) / units->year, /* softness parameter (Pa^{-n}/s) */
368: B = PetscPowReal(A,-1./n); /* hardness parameter */
369: thi->viscosity.Bd2 = B/2;
370: thi->viscosity.exponent = (p-2)/2;
371: thi->viscosity.eps = 0.5*PetscSqr(1e-5 / units->year);
372: }
373: Bd2 = thi->viscosity.Bd2;
374: exponent = thi->viscosity.exponent;
375: eps = thi->viscosity.eps;
376: *eta = Bd2 * PetscPowReal(eps + gam,exponent);
377: *deta = exponent * (*eta) / (eps + gam);
378: }
380: static void RangeUpdate(PetscReal *min,PetscReal *max,PetscReal x)
381: {
382: if (x < *min) *min = x;
383: if (x > *max) *max = x;
384: }
386: static void PRangeClear(PRange *p)
387: {
388: p->cmin = p->min = 1e100;
389: p->cmax = p->max = -1e100;
390: }
392: static PetscErrorCode PRangeMinMax(PRange *p,PetscReal min,PetscReal max)
393: {
396: p->cmin = min;
397: p->cmax = max;
398: if (min < p->min) p->min = min;
399: if (max > p->max) p->max = max;
400: return(0);
401: }
403: static PetscErrorCode THIDestroy(THI *thi)
404: {
408: if (!*thi) return(0);
409: if (--((PetscObject)(*thi))->refct > 0) {*thi = 0; return(0);}
410: PetscFree((*thi)->units);
411: PetscFree((*thi)->mattype);
412: PetscHeaderDestroy(thi);
413: return(0);
414: }
416: static PetscErrorCode THICreate(MPI_Comm comm,THI *inthi)
417: {
418: static PetscBool registered = PETSC_FALSE;
419: THI thi;
420: Units units;
421: PetscErrorCode ierr;
424: *inthi = 0;
425: if (!registered) {
426: PetscClassIdRegister("Toy Hydrostatic Ice",&THI_CLASSID);
427: registered = PETSC_TRUE;
428: }
429: PetscHeaderCreate(thi,THI_CLASSID,"THI","Toy Hydrostatic Ice","",comm,THIDestroy,0);
431: PetscNew(&thi->units);
432: units = thi->units;
433: units->meter = 1e-2;
434: units->second = 1e-7;
435: units->kilogram = 1e-12;
437: PetscOptionsBegin(comm,NULL,"Scaled units options","");
438: {
439: PetscOptionsReal("-units_meter","1 meter in scaled length units","",units->meter,&units->meter,NULL);
440: PetscOptionsReal("-units_second","1 second in scaled time units","",units->second,&units->second,NULL);
441: PetscOptionsReal("-units_kilogram","1 kilogram in scaled mass units","",units->kilogram,&units->kilogram,NULL);
442: }
443: PetscOptionsEnd();
444: units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second));
445: units->year = 31556926. * units->second; /* seconds per year */
447: thi->Lx = 10.e3;
448: thi->Ly = 10.e3;
449: thi->Lz = 1000;
450: thi->dirichlet_scale = 1;
451: thi->verbose = PETSC_FALSE;
453: PetscOptionsBegin(comm,NULL,"Toy Hydrostatic Ice options","");
454: {
455: QuadratureType quad = QUAD_GAUSS;
456: char homexp[] = "A";
457: char mtype[256] = MATSBAIJ;
458: PetscReal L,m = 1.0;
459: PetscBool flg;
460: L = thi->Lx;
461: PetscOptionsReal("-thi_L","Domain size (m)","",L,&L,&flg);
462: if (flg) thi->Lx = thi->Ly = L;
463: PetscOptionsReal("-thi_Lx","X Domain size (m)","",thi->Lx,&thi->Lx,NULL);
464: PetscOptionsReal("-thi_Ly","Y Domain size (m)","",thi->Ly,&thi->Ly,NULL);
465: PetscOptionsReal("-thi_Lz","Z Domain size (m)","",thi->Lz,&thi->Lz,NULL);
466: PetscOptionsString("-thi_hom","ISMIP-HOM experiment (A or C)","",homexp,homexp,sizeof(homexp),NULL);
467: switch (homexp[0] = toupper(homexp[0])) {
468: case 'A':
469: thi->initialize = THIInitialize_HOM_A;
470: thi->no_slip = PETSC_TRUE;
471: thi->alpha = 0.5;
472: break;
473: case 'C':
474: thi->initialize = THIInitialize_HOM_C;
475: thi->no_slip = PETSC_FALSE;
476: thi->alpha = 0.1;
477: break;
478: case 'X':
479: thi->initialize = THIInitialize_HOM_X;
480: thi->no_slip = PETSC_FALSE;
481: thi->alpha = 0.3;
482: break;
483: case 'Y':
484: thi->initialize = THIInitialize_HOM_Y;
485: thi->no_slip = PETSC_FALSE;
486: thi->alpha = 0.5;
487: break;
488: case 'Z':
489: thi->initialize = THIInitialize_HOM_Z;
490: thi->no_slip = PETSC_FALSE;
491: thi->alpha = 0.5;
492: break;
493: default:
494: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"HOM experiment '%c' not implemented",homexp[0]);
495: }
496: PetscOptionsEnum("-thi_quadrature","Quadrature to use for 3D elements","",QuadratureTypes,(PetscEnum)quad,(PetscEnum*)&quad,NULL);
497: switch (quad) {
498: case QUAD_GAUSS:
499: HexQInterp = HexQInterp_Gauss;
500: HexQDeriv = HexQDeriv_Gauss;
501: break;
502: case QUAD_LOBATTO:
503: HexQInterp = HexQInterp_Lobatto;
504: HexQDeriv = HexQDeriv_Lobatto;
505: break;
506: }
507: PetscOptionsReal("-thi_alpha","Bed angle (degrees)","",thi->alpha,&thi->alpha,NULL);
509: thi->friction.refvel = 100.;
510: thi->friction.epsvel = 1.;
512: PetscOptionsReal("-thi_friction_refvel","Reference velocity for sliding","",thi->friction.refvel,&thi->friction.refvel,NULL);
513: PetscOptionsReal("-thi_friction_epsvel","Regularization velocity for sliding","",thi->friction.epsvel,&thi->friction.epsvel,NULL);
514: PetscOptionsReal("-thi_friction_m","Friction exponent, 0=Coulomb, 1=Navier","",m,&m,NULL);
516: thi->friction.exponent = (m-1)/2;
518: PetscOptionsReal("-thi_dirichlet_scale","Scale Dirichlet boundary conditions by this factor","",thi->dirichlet_scale,&thi->dirichlet_scale,NULL);
519: 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);
520: PetscOptionsBool("-thi_coarse2d","Use a 2D coarse space corresponding to SSA","",thi->coarse2d,&thi->coarse2d,NULL);
521: PetscOptionsBool("-thi_tridiagonal","Assemble a tridiagonal system (column coupling only) on the finest level","",thi->tridiagonal,&thi->tridiagonal,NULL);
522: PetscOptionsFList("-thi_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)mtype,sizeof(mtype),NULL);
523: PetscStrallocpy(mtype,(char**)&thi->mattype);
524: PetscOptionsBool("-thi_verbose","Enable verbose output (like matrix sizes and statistics)","",thi->verbose,&thi->verbose,NULL);
525: }
526: PetscOptionsEnd();
528: /* dimensionalize */
529: thi->Lx *= units->meter;
530: thi->Ly *= units->meter;
531: thi->Lz *= units->meter;
532: thi->alpha *= PETSC_PI / 180;
534: PRangeClear(&thi->eta);
535: PRangeClear(&thi->beta2);
537: {
538: PetscReal u = 1000*units->meter/(3e7*units->second),
539: gradu = u / (100*units->meter),eta,deta,
540: rho = 910 * units->kilogram/PetscPowReal(units->meter,3),
541: grav = 9.81 * units->meter/PetscSqr(units->second),
542: driving = rho * grav * PetscSinReal(thi->alpha) * 1000*units->meter;
543: THIViscosity(thi,0.5*gradu*gradu,&eta,&deta);
544: thi->rhog = rho * grav;
545: if (thi->verbose) {
546: 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);
547: 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);
548: 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));
549: THIViscosity(thi,0.5*PetscSqr(1e-3*gradu),&eta,&deta);
550: 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));
551: }
552: }
554: *inthi = thi;
555: return(0);
556: }
558: static PetscErrorCode THIInitializePrm(THI thi,DM da2prm,Vec prm)
559: {
560: PrmNode **p;
561: PetscInt i,j,xs,xm,ys,ym,mx,my;
565: DMDAGetGhostCorners(da2prm,&ys,&xs,0,&ym,&xm,0);
566: DMDAGetInfo(da2prm,0, &my,&mx,0, 0,0,0, 0,0,0,0,0,0);
567: DMDAVecGetArray(da2prm,prm,&p);
568: for (i=xs; i<xs+xm; i++) {
569: for (j=ys; j<ys+ym; j++) {
570: PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my;
571: thi->initialize(thi,xx,yy,&p[i][j]);
572: }
573: }
574: DMDAVecRestoreArray(da2prm,prm,&p);
575: return(0);
576: }
578: static PetscErrorCode THISetUpDM(THI thi,DM dm)
579: {
580: PetscErrorCode ierr;
581: PetscInt refinelevel,coarsenlevel,level,dim,Mx,My,Mz,mx,my,s;
582: DMDAStencilType st;
583: DM da2prm;
584: Vec X;
587: DMDAGetInfo(dm,&dim, &Mz,&My,&Mx, 0,&my,&mx, 0,&s,0,0,0,&st);
588: if (dim == 2) {
589: DMDAGetInfo(dm,&dim, &My,&Mx,0, &my,&mx,0, 0,&s,0,0,0,&st);
590: }
591: DMGetRefineLevel(dm,&refinelevel);
592: DMGetCoarsenLevel(dm,&coarsenlevel);
593: level = refinelevel - coarsenlevel;
594: DMDACreate2d(PetscObjectComm((PetscObject)thi),DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,My,Mx,my,mx,sizeof(PrmNode)/sizeof(PetscScalar),s,0,0,&da2prm);
595: DMSetUp(da2prm);
596: DMCreateLocalVector(da2prm,&X);
597: {
598: PetscReal Lx = thi->Lx / thi->units->meter,Ly = thi->Ly / thi->units->meter,Lz = thi->Lz / thi->units->meter;
599: if (dim == 2) {
600: 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));
601: } else {
602: 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)));
603: }
604: }
605: THIInitializePrm(thi,da2prm,X);
606: if (thi->tridiagonal) { /* Reset coarse Jacobian evaluation */
607: DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_3D_Full,thi);
608: }
609: if (thi->coarse2d) {
610: DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_2D,thi);
611: }
612: PetscObjectCompose((PetscObject)dm,"DMDA2Prm",(PetscObject)da2prm);
613: PetscObjectCompose((PetscObject)dm,"DMDA2Prm_Vec",(PetscObject)X);
614: DMDestroy(&da2prm);
615: VecDestroy(&X);
616: return(0);
617: }
619: static PetscErrorCode DMCoarsenHook_THI(DM dmf,DM dmc,void *ctx)
620: {
621: THI thi = (THI)ctx;
623: PetscInt rlevel,clevel;
626: THISetUpDM(thi,dmc);
627: DMGetRefineLevel(dmc,&rlevel);
628: DMGetCoarsenLevel(dmc,&clevel);
629: if (rlevel-clevel == 0) {DMSetMatType(dmc,MATAIJ);}
630: DMCoarsenHookAdd(dmc,DMCoarsenHook_THI,NULL,thi);
631: return(0);
632: }
634: static PetscErrorCode DMRefineHook_THI(DM dmc,DM dmf,void *ctx)
635: {
636: THI thi = (THI)ctx;
640: THISetUpDM(thi,dmf);
641: DMSetMatType(dmf,thi->mattype);
642: DMRefineHookAdd(dmf,DMRefineHook_THI,NULL,thi);
643: /* With grid sequencing, a formerly-refined DM will later be coarsened by PCSetUp_MG */
644: DMCoarsenHookAdd(dmf,DMCoarsenHook_THI,NULL,thi);
645: return(0);
646: }
648: static PetscErrorCode THIDAGetPrm(DM da,PrmNode ***prm)
649: {
651: DM da2prm;
652: Vec X;
655: PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
656: if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
657: PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
658: if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
659: DMDAVecGetArray(da2prm,X,prm);
660: return(0);
661: }
663: static PetscErrorCode THIDARestorePrm(DM da,PrmNode ***prm)
664: {
666: DM da2prm;
667: Vec X;
670: PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
671: if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
672: PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
673: if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
674: DMDAVecRestoreArray(da2prm,X,prm);
675: return(0);
676: }
678: static PetscErrorCode THIInitial(SNES snes,Vec X,void *ctx)
679: {
680: THI thi;
681: PetscInt i,j,k,xs,xm,ys,ym,zs,zm,mx,my;
682: PetscReal hx,hy;
683: PrmNode **prm;
684: Node ***x;
686: DM da;
689: SNESGetDM(snes,&da);
690: DMGetApplicationContext(da,&thi);
691: DMDAGetInfo(da,0, 0,&my,&mx, 0,0,0, 0,0,0,0,0,0);
692: DMDAGetCorners(da,&zs,&ys,&xs,&zm,&ym,&xm);
693: DMDAVecGetArray(da,X,&x);
694: THIDAGetPrm(da,&prm);
695: hx = thi->Lx / mx;
696: hy = thi->Ly / my;
697: for (i=xs; i<xs+xm; i++) {
698: for (j=ys; j<ys+ym; j++) {
699: for (k=zs; k<zs+zm; k++) {
700: const PetscScalar zm1 = zm-1,
701: 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),
702: 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);
703: x[i][j][k].u = 0. * drivingx * prm[i][j].h*(PetscScalar)k/zm1;
704: x[i][j][k].v = 0. * drivingy * prm[i][j].h*(PetscScalar)k/zm1;
705: }
706: }
707: }
708: DMDAVecRestoreArray(da,X,&x);
709: THIDARestorePrm(da,&prm);
710: return(0);
711: }
713: 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)
714: {
715: PetscInt l,ll;
716: PetscScalar gam;
718: du[0] = du[1] = du[2] = 0;
719: dv[0] = dv[1] = dv[2] = 0;
720: *u = 0;
721: *v = 0;
722: for (l=0; l<8; l++) {
723: *u += phi[l] * n[l].u;
724: *v += phi[l] * n[l].v;
725: for (ll=0; ll<3; ll++) {
726: du[ll] += dphi[l][ll] * n[l].u;
727: dv[ll] += dphi[l][ll] * n[l].v;
728: }
729: }
730: 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]);
731: THIViscosity(thi,PetscRealPart(gam),eta,deta);
732: }
734: 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)
735: {
736: PetscInt l,ll;
737: PetscScalar gam;
739: du[0] = du[1] = 0;
740: dv[0] = dv[1] = 0;
741: *u = 0;
742: *v = 0;
743: for (l=0; l<4; l++) {
744: *u += phi[l] * n[l].u;
745: *v += phi[l] * n[l].v;
746: for (ll=0; ll<2; ll++) {
747: du[ll] += dphi[l][ll] * n[l].u;
748: dv[ll] += dphi[l][ll] * n[l].v;
749: }
750: }
751: gam = PetscSqr(du[0]) + PetscSqr(dv[1]) + du[0]*dv[1] + 0.25*PetscSqr(du[1]+dv[0]);
752: THIViscosity(thi,PetscRealPart(gam),eta,deta);
753: }
755: static PetscErrorCode THIFunctionLocal(DMDALocalInfo *info,Node ***x,Node ***f,THI thi)
756: {
757: PetscInt xs,ys,xm,ym,zm,i,j,k,q,l;
758: PetscReal hx,hy,etamin,etamax,beta2min,beta2max;
759: PrmNode **prm;
763: xs = info->zs;
764: ys = info->ys;
765: xm = info->zm;
766: ym = info->ym;
767: zm = info->xm;
768: hx = thi->Lx / info->mz;
769: hy = thi->Ly / info->my;
771: etamin = 1e100;
772: etamax = 0;
773: beta2min = 1e100;
774: beta2max = 0;
776: THIDAGetPrm(info->da,&prm);
778: for (i=xs; i<xs+xm; i++) {
779: for (j=ys; j<ys+ym; j++) {
780: PrmNode pn[4];
781: QuadExtract(prm,i,j,pn);
782: for (k=0; k<zm-1; k++) {
783: PetscInt ls = 0;
784: Node n[8],*fn[8];
785: PetscReal zn[8],etabase = 0;
786: PrmHexGetZ(pn,k,zm,zn);
787: HexExtract(x,i,j,k,n);
788: HexExtractRef(f,i,j,k,fn);
789: if (thi->no_slip && k == 0) {
790: for (l=0; l<4; l++) n[l].u = n[l].v = 0;
791: /* The first 4 basis functions lie on the bottom layer, so their contribution is exactly 0, hence we can skip them */
792: ls = 4;
793: }
794: for (q=0; q<8; q++) {
795: PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta;
796: PetscScalar du[3],dv[3],u,v;
797: HexGrad(HexQDeriv[q],zn,dz);
798: HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
799: PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
800: jw /= thi->rhog; /* scales residuals to be O(1) */
801: if (q == 0) etabase = eta;
802: RangeUpdate(&etamin,&etamax,eta);
803: for (l=ls; l<8; l++) { /* test functions */
804: const PetscReal ds[2] = {-PetscSinReal(thi->alpha),0};
805: const PetscReal pp = phi[l],*dp = dphi[l];
806: 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];
807: 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];
808: }
809: }
810: if (k == 0) { /* we are on a bottom face */
811: if (thi->no_slip) {
812: /* Note: Non-Galerkin coarse grid operators are very sensitive to the scaling of Dirichlet boundary
813: * conditions. After shenanigans above, etabase contains the effective viscosity at the closest quadrature
814: * point to the bed. We want the diagonal entry in the Dirichlet condition to have similar magnitude to the
815: * diagonal entry corresponding to the adjacent node. The fundamental scaling of the viscous part is in
816: * diagu, diagv below. This scaling is easy to recognize by considering the finite difference operator after
817: * scaling by element size. The no-slip Dirichlet condition is scaled by this factor, and also in the
818: * assembled matrix (see the similar block in THIJacobianLocal).
819: *
820: * Note that the residual at this Dirichlet node is linear in the state at this node, but also depends
821: * (nonlinearly in general) on the neighboring interior nodes through the local viscosity. This will make
822: * a matrix-free Jacobian have extra entries in the corresponding row. We assemble only the diagonal part,
823: * so the solution will exactly satisfy the boundary condition after the first linear iteration.
824: */
825: const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1.);
826: 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);
827: fn[0]->u = thi->dirichlet_scale*diagu*x[i][j][k].u;
828: fn[0]->v = thi->dirichlet_scale*diagv*x[i][j][k].v;
829: } else { /* Integrate over bottom face to apply boundary condition */
830: for (q=0; q<4; q++) {
831: const PetscReal jw = 0.25*hx*hy/thi->rhog,*phi = QuadQInterp[q];
832: PetscScalar u =0,v=0,rbeta2=0;
833: PetscReal beta2,dbeta2;
834: for (l=0; l<4; l++) {
835: u += phi[l]*n[l].u;
836: v += phi[l]*n[l].v;
837: rbeta2 += phi[l]*pn[l].beta2;
838: }
839: THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
840: RangeUpdate(&beta2min,&beta2max,beta2);
841: for (l=0; l<4; l++) {
842: const PetscReal pp = phi[l];
843: fn[ls+l]->u += pp*jw*beta2*u;
844: fn[ls+l]->v += pp*jw*beta2*v;
845: }
846: }
847: }
848: }
849: }
850: }
851: }
853: THIDARestorePrm(info->da,&prm);
855: PRangeMinMax(&thi->eta,etamin,etamax);
856: PRangeMinMax(&thi->beta2,beta2min,beta2max);
857: return(0);
858: }
860: static PetscErrorCode THIMatrixStatistics(THI thi,Mat B,PetscViewer viewer)
861: {
863: PetscReal nrm;
864: PetscInt m;
865: PetscMPIInt rank;
868: MatNorm(B,NORM_FROBENIUS,&nrm);
869: MatGetSize(B,&m,0);
870: MPI_Comm_rank(PetscObjectComm((PetscObject)B),&rank);
871: if (!rank) {
872: PetscScalar val0,val2;
873: MatGetValue(B,0,0,&val0);
874: MatGetValue(B,2,2,&val2);
875: 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);
876: }
877: return(0);
878: }
880: static PetscErrorCode THISurfaceStatistics(DM da,Vec X,PetscReal *min,PetscReal *max,PetscReal *mean)
881: {
883: Node ***x;
884: PetscInt i,j,xs,ys,zs,xm,ym,zm,mx,my,mz;
885: PetscReal umin = 1e100,umax=-1e100;
886: PetscScalar usum = 0.0,gusum;
889: *min = *max = *mean = 0;
890: DMDAGetInfo(da,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
891: DMDAGetCorners(da,&zs,&ys,&xs,&zm,&ym,&xm);
892: if (zs != 0 || zm != mz) SETERRQ(PETSC_COMM_SELF,1,"Unexpected decomposition");
893: DMDAVecGetArray(da,X,&x);
894: for (i=xs; i<xs+xm; i++) {
895: for (j=ys; j<ys+ym; j++) {
896: PetscReal u = PetscRealPart(x[i][j][zm-1].u);
897: RangeUpdate(&umin,&umax,u);
898: usum += u;
899: }
900: }
901: DMDAVecRestoreArray(da,X,&x);
902: MPI_Allreduce(&umin,min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));
903: MPI_Allreduce(&umax,max,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));
904: MPI_Allreduce(&usum,&gusum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da));
905: *mean = PetscRealPart(gusum) / (mx*my);
906: return(0);
907: }
909: static PetscErrorCode THISolveStatistics(THI thi,SNES snes,PetscInt coarsened,const char name[])
910: {
911: MPI_Comm comm;
912: Vec X;
913: DM dm;
917: PetscObjectGetComm((PetscObject)thi,&comm);
918: SNESGetSolution(snes,&X);
919: SNESGetDM(snes,&dm);
920: PetscPrintf(comm,"Solution statistics after solve: %s\n",name);
921: {
922: PetscInt its,lits;
923: SNESConvergedReason reason;
924: SNESGetIterationNumber(snes,&its);
925: SNESGetConvergedReason(snes,&reason);
926: SNESGetLinearSolveIterations(snes,&lits);
927: PetscPrintf(comm,"%s: Number of SNES iterations = %D, total linear iterations = %D\n",SNESConvergedReasons[reason],its,lits);
928: }
929: {
930: PetscReal nrm2,tmin[3]={1e100,1e100,1e100},tmax[3]={-1e100,-1e100,-1e100},min[3],max[3];
931: PetscInt i,j,m;
932: const PetscScalar *x;
933: VecNorm(X,NORM_2,&nrm2);
934: VecGetLocalSize(X,&m);
935: VecGetArrayRead(X,&x);
936: for (i=0; i<m; i+=2) {
937: PetscReal u = PetscRealPart(x[i]),v = PetscRealPart(x[i+1]),c = PetscSqrtReal(u*u+v*v);
938: tmin[0] = PetscMin(u,tmin[0]);
939: tmin[1] = PetscMin(v,tmin[1]);
940: tmin[2] = PetscMin(c,tmin[2]);
941: tmax[0] = PetscMax(u,tmax[0]);
942: tmax[1] = PetscMax(v,tmax[1]);
943: tmax[2] = PetscMax(c,tmax[2]);
944: }
945: VecRestoreArrayRead(X,&x);
946: MPI_Allreduce(tmin,min,3,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)thi));
947: MPI_Allreduce(tmax,max,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)thi));
948: /* Dimensionalize to meters/year */
949: nrm2 *= thi->units->year / thi->units->meter;
950: for (j=0; j<3; j++) {
951: min[j] *= thi->units->year / thi->units->meter;
952: max[j] *= thi->units->year / thi->units->meter;
953: }
954: if (min[0] == 0.0) min[0] = 0.0;
955: 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]);
956: {
957: PetscReal umin,umax,umean;
958: THISurfaceStatistics(dm,X,&umin,&umax,&umean);
959: umin *= thi->units->year / thi->units->meter;
960: umax *= thi->units->year / thi->units->meter;
961: umean *= thi->units->year / thi->units->meter;
962: PetscPrintf(comm,"Surface statistics: u in [%12.6e, %12.6e] mean %12.6e\n",(double)umin,(double)umax,(double)umean);
963: }
964: /* These values stay nondimensional */
965: 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);
966: 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);
967: }
968: return(0);
969: }
971: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info,Node **x,Mat J,Mat B,THI thi)
972: {
973: PetscInt xs,ys,xm,ym,i,j,q,l,ll;
974: PetscReal hx,hy;
975: PrmNode **prm;
979: xs = info->ys;
980: ys = info->xs;
981: xm = info->ym;
982: ym = info->xm;
983: hx = thi->Lx / info->my;
984: hy = thi->Ly / info->mx;
986: MatZeroEntries(B);
987: THIDAGetPrm(info->da,&prm);
989: for (i=xs; i<xs+xm; i++) {
990: for (j=ys; j<ys+ym; j++) {
991: Node n[4];
992: PrmNode pn[4];
993: PetscScalar Ke[4*2][4*2];
994: QuadExtract(prm,i,j,pn);
995: QuadExtract(x,i,j,n);
996: PetscMemzero(Ke,sizeof(Ke));
997: for (q=0; q<4; q++) {
998: PetscReal phi[4],dphi[4][2],jw,eta,deta,beta2,dbeta2;
999: PetscScalar u,v,du[2],dv[2],h = 0,rbeta2 = 0;
1000: for (l=0; l<4; l++) {
1001: phi[l] = QuadQInterp[q][l];
1002: dphi[l][0] = QuadQDeriv[q][l][0]*2./hx;
1003: dphi[l][1] = QuadQDeriv[q][l][1]*2./hy;
1004: h += phi[l] * pn[l].h;
1005: rbeta2 += phi[l] * pn[l].beta2;
1006: }
1007: jw = 0.25*hx*hy / thi->rhog; /* rhog is only scaling */
1008: PointwiseNonlinearity2D(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
1009: THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
1010: for (l=0; l<4; l++) {
1011: const PetscReal pp = phi[l],*dp = dphi[l];
1012: for (ll=0; ll<4; ll++) {
1013: const PetscReal ppl = phi[ll],*dpl = dphi[ll];
1014: PetscScalar dgdu,dgdv;
1015: dgdu = 2.*du[0]*dpl[0] + dv[1]*dpl[0] + 0.5*(du[1]+dv[0])*dpl[1];
1016: dgdv = 2.*dv[1]*dpl[1] + du[0]*dpl[1] + 0.5*(du[1]+dv[0])*dpl[0];
1017: /* Picard part */
1018: 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;
1019: Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0];
1020: Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1];
1021: 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;
1022: /* extra Newton terms */
1023: 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;
1024: 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;
1025: 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;
1026: 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;
1027: }
1028: }
1029: }
1030: {
1031: const MatStencil rc[4] = {{0,i,j,0},{0,i+1,j,0},{0,i+1,j+1,0},{0,i,j+1,0}};
1032: MatSetValuesBlockedStencil(B,4,rc,4,rc,&Ke[0][0],ADD_VALUES);
1033: }
1034: }
1035: }
1036: THIDARestorePrm(info->da,&prm);
1038: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1039: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1040: MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1041: if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1042: return(0);
1043: }
1045: static PetscErrorCode THIJacobianLocal_3D(DMDALocalInfo *info,Node ***x,Mat B,THI thi,THIAssemblyMode amode)
1046: {
1047: PetscInt xs,ys,xm,ym,zm,i,j,k,q,l,ll;
1048: PetscReal hx,hy;
1049: PrmNode **prm;
1053: xs = info->zs;
1054: ys = info->ys;
1055: xm = info->zm;
1056: ym = info->ym;
1057: zm = info->xm;
1058: hx = thi->Lx / info->mz;
1059: hy = thi->Ly / info->my;
1061: MatZeroEntries(B);
1062: MatSetOption(B,MAT_SUBSET_OFF_PROC_ENTRIES,PETSC_TRUE);
1063: THIDAGetPrm(info->da,&prm);
1065: for (i=xs; i<xs+xm; i++) {
1066: for (j=ys; j<ys+ym; j++) {
1067: PrmNode pn[4];
1068: QuadExtract(prm,i,j,pn);
1069: for (k=0; k<zm-1; k++) {
1070: Node n[8];
1071: PetscReal zn[8],etabase = 0;
1072: PetscScalar Ke[8*2][8*2];
1073: PetscInt ls = 0;
1075: PrmHexGetZ(pn,k,zm,zn);
1076: HexExtract(x,i,j,k,n);
1077: PetscMemzero(Ke,sizeof(Ke));
1078: if (thi->no_slip && k == 0) {
1079: for (l=0; l<4; l++) n[l].u = n[l].v = 0;
1080: ls = 4;
1081: }
1082: for (q=0; q<8; q++) {
1083: PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta;
1084: PetscScalar du[3],dv[3],u,v;
1085: HexGrad(HexQDeriv[q],zn,dz);
1086: HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
1087: PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
1088: jw /= thi->rhog; /* residuals are scaled by this factor */
1089: if (q == 0) etabase = eta;
1090: for (l=ls; l<8; l++) { /* test functions */
1091: const PetscReal *PETSC_RESTRICT dp = dphi[l];
1092: #if USE_SSE2_KERNELS
1093: /* gcc (up to my 4.5 snapshot) is really bad at hoisting intrinsics so we do it manually */
1094: __m128d
1095: p4 = _mm_set1_pd(4),p2 = _mm_set1_pd(2),p05 = _mm_set1_pd(0.5),
1096: p42 = _mm_setr_pd(4,2),p24 = _mm_shuffle_pd(p42,p42,_MM_SHUFFLE2(0,1)),
1097: du0 = _mm_set1_pd(du[0]),du1 = _mm_set1_pd(du[1]),du2 = _mm_set1_pd(du[2]),
1098: dv0 = _mm_set1_pd(dv[0]),dv1 = _mm_set1_pd(dv[1]),dv2 = _mm_set1_pd(dv[2]),
1099: jweta = _mm_set1_pd(jw*eta),jwdeta = _mm_set1_pd(jw*deta),
1100: dp0 = _mm_set1_pd(dp[0]),dp1 = _mm_set1_pd(dp[1]),dp2 = _mm_set1_pd(dp[2]),
1101: dp0jweta = _mm_mul_pd(dp0,jweta),dp1jweta = _mm_mul_pd(dp1,jweta),dp2jweta = _mm_mul_pd(dp2,jweta),
1102: p4du0p2dv1 = _mm_add_pd(_mm_mul_pd(p4,du0),_mm_mul_pd(p2,dv1)), /* 4 du0 + 2 dv1 */
1103: p4dv1p2du0 = _mm_add_pd(_mm_mul_pd(p4,dv1),_mm_mul_pd(p2,du0)), /* 4 dv1 + 2 du0 */
1104: pdu2dv2 = _mm_unpacklo_pd(du2,dv2), /* [du2, dv2] */
1105: du1pdv0 = _mm_add_pd(du1,dv0), /* du1 + dv0 */
1106: t1 = _mm_mul_pd(dp0,p4du0p2dv1), /* dp0 (4 du0 + 2 dv1) */
1107: t2 = _mm_mul_pd(dp1,p4dv1p2du0); /* dp1 (4 dv1 + 2 du0) */
1109: #endif
1110: #if defined COMPUTE_LOWER_TRIANGULAR /* The element matrices are always symmetric so computing the lower-triangular part is not necessary */
1111: for (ll=ls; ll<8; ll++) { /* trial functions */
1112: #else
1113: for (ll=l; ll<8; ll++) {
1114: #endif
1115: const PetscReal *PETSC_RESTRICT dpl = dphi[ll];
1116: if (amode == THIASSEMBLY_TRIDIAGONAL && (l-ll)%4) continue; /* these entries would not be inserted */
1117: #if !USE_SSE2_KERNELS
1118: /* The analytic Jacobian in nice, easy-to-read form */
1119: {
1120: PetscScalar dgdu,dgdv;
1121: dgdu = 2.*du[0]*dpl[0] + dv[1]*dpl[0] + 0.5*(du[1]+dv[0])*dpl[1] + 0.5*du[2]*dpl[2];
1122: dgdv = 2.*dv[1]*dpl[1] + du[0]*dpl[1] + 0.5*(du[1]+dv[0])*dpl[0] + 0.5*dv[2]*dpl[2];
1123: /* Picard part */
1124: 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];
1125: Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0];
1126: Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1];
1127: 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];
1128: /* extra Newton terms */
1129: 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];
1130: 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];
1131: 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];
1132: 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];
1133: }
1134: #else
1135: /* This SSE2 code is an exact replica of above, but uses explicit packed instructions for some speed
1136: * benefit. On my hardware, these intrinsics are almost twice as fast as above, reducing total assembly cost
1137: * by 25 to 30 percent. */
1138: {
1139: __m128d
1140: keu = _mm_loadu_pd(&Ke[l*2+0][ll*2+0]),
1141: kev = _mm_loadu_pd(&Ke[l*2+1][ll*2+0]),
1142: dpl01 = _mm_loadu_pd(&dpl[0]),dpl10 = _mm_shuffle_pd(dpl01,dpl01,_MM_SHUFFLE2(0,1)),dpl2 = _mm_set_sd(dpl[2]),
1143: t0,t3,pdgduv;
1144: keu = _mm_add_pd(keu,_mm_add_pd(_mm_mul_pd(_mm_mul_pd(dp0jweta,p42),dpl01),
1145: _mm_add_pd(_mm_mul_pd(dp1jweta,dpl10),
1146: _mm_mul_pd(dp2jweta,dpl2))));
1147: kev = _mm_add_pd(kev,_mm_add_pd(_mm_mul_pd(_mm_mul_pd(dp1jweta,p24),dpl01),
1148: _mm_add_pd(_mm_mul_pd(dp0jweta,dpl10),
1149: _mm_mul_pd(dp2jweta,_mm_shuffle_pd(dpl2,dpl2,_MM_SHUFFLE2(0,1))))));
1150: pdgduv = _mm_mul_pd(p05,_mm_add_pd(_mm_add_pd(_mm_mul_pd(p42,_mm_mul_pd(du0,dpl01)),
1151: _mm_mul_pd(p24,_mm_mul_pd(dv1,dpl01))),
1152: _mm_add_pd(_mm_mul_pd(du1pdv0,dpl10),
1153: _mm_mul_pd(pdu2dv2,_mm_set1_pd(dpl[2]))))); /* [dgdu, dgdv] */
1154: t0 = _mm_mul_pd(jwdeta,pdgduv); /* jw deta [dgdu, dgdv] */
1155: t3 = _mm_mul_pd(t0,du1pdv0); /* t0 (du1 + dv0) */
1156: _mm_storeu_pd(&Ke[l*2+0][ll*2+0],_mm_add_pd(keu,_mm_add_pd(_mm_mul_pd(t1,t0),
1157: _mm_add_pd(_mm_mul_pd(dp1,t3),
1158: _mm_mul_pd(t0,_mm_mul_pd(dp2,du2))))));
1159: _mm_storeu_pd(&Ke[l*2+1][ll*2+0],_mm_add_pd(kev,_mm_add_pd(_mm_mul_pd(t2,t0),
1160: _mm_add_pd(_mm_mul_pd(dp0,t3),
1161: _mm_mul_pd(t0,_mm_mul_pd(dp2,dv2))))));
1162: }
1163: #endif
1164: }
1165: }
1166: }
1167: if (k == 0) { /* on a bottom face */
1168: if (thi->no_slip) {
1169: const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1);
1170: 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);
1171: Ke[0][0] = thi->dirichlet_scale*diagu;
1172: Ke[1][1] = thi->dirichlet_scale*diagv;
1173: } else {
1174: for (q=0; q<4; q++) {
1175: const PetscReal jw = 0.25*hx*hy/thi->rhog,*phi = QuadQInterp[q];
1176: PetscScalar u =0,v=0,rbeta2=0;
1177: PetscReal beta2,dbeta2;
1178: for (l=0; l<4; l++) {
1179: u += phi[l]*n[l].u;
1180: v += phi[l]*n[l].v;
1181: rbeta2 += phi[l]*pn[l].beta2;
1182: }
1183: THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
1184: for (l=0; l<4; l++) {
1185: const PetscReal pp = phi[l];
1186: for (ll=0; ll<4; ll++) {
1187: const PetscReal ppl = phi[ll];
1188: Ke[l*2+0][ll*2+0] += pp*jw*beta2*ppl + pp*jw*dbeta2*u*u*ppl;
1189: Ke[l*2+0][ll*2+1] += pp*jw*dbeta2*u*v*ppl;
1190: Ke[l*2+1][ll*2+0] += pp*jw*dbeta2*v*u*ppl;
1191: Ke[l*2+1][ll*2+1] += pp*jw*beta2*ppl + pp*jw*dbeta2*v*v*ppl;
1192: }
1193: }
1194: }
1195: }
1196: }
1197: {
1198: 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}};
1199: if (amode == THIASSEMBLY_TRIDIAGONAL) {
1200: for (l=0; l<4; l++) { /* Copy out each of the blocks, discarding horizontal coupling */
1201: const PetscInt l4 = l+4;
1202: const MatStencil rcl[2] = {{rc[l].k,rc[l].j,rc[l].i,0},{rc[l4].k,rc[l4].j,rc[l4].i,0}};
1203: #if defined COMPUTE_LOWER_TRIANGULAR
1204: 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]},
1205: {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]},
1206: {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]},
1207: {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]}};
1208: #else
1209: /* Same as above except for the lower-left block */
1210: 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]},
1211: {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]},
1212: {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]},
1213: {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]}};
1214: #endif
1215: MatSetValuesBlockedStencil(B,2,rcl,2,rcl,&Kel[0][0],ADD_VALUES);
1216: }
1217: } else {
1218: #if !defined COMPUTE_LOWER_TRIANGULAR /* fill in lower-triangular part, this is really cheap compared to computing the entries */
1219: for (l=0; l<8; l++) {
1220: for (ll=l+1; ll<8; ll++) {
1221: Ke[ll*2+0][l*2+0] = Ke[l*2+0][ll*2+0];
1222: Ke[ll*2+1][l*2+0] = Ke[l*2+0][ll*2+1];
1223: Ke[ll*2+0][l*2+1] = Ke[l*2+1][ll*2+0];
1224: Ke[ll*2+1][l*2+1] = Ke[l*2+1][ll*2+1];
1225: }
1226: }
1227: #endif
1228: MatSetValuesBlockedStencil(B,8,rc,8,rc,&Ke[0][0],ADD_VALUES);
1229: }
1230: }
1231: }
1232: }
1233: }
1234: THIDARestorePrm(info->da,&prm);
1236: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1237: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1238: MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1239: if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1240: return(0);
1241: }
1243: static PetscErrorCode THIJacobianLocal_3D_Full(DMDALocalInfo *info,Node ***x,Mat A,Mat B,THI thi)
1244: {
1248: THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_FULL);
1249: return(0);
1250: }
1252: static PetscErrorCode THIJacobianLocal_3D_Tridiagonal(DMDALocalInfo *info,Node ***x,Mat A,Mat B,THI thi)
1253: {
1257: THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_TRIDIAGONAL);
1258: return(0);
1259: }
1261: static PetscErrorCode DMRefineHierarchy_THI(DM dac0,PetscInt nlevels,DM hierarchy[])
1262: {
1263: PetscErrorCode ierr;
1264: THI thi;
1265: PetscInt dim,M,N,m,n,s,dof;
1266: DM dac,daf;
1267: DMDAStencilType st;
1268: DM_DA *ddf,*ddc;
1271: PetscObjectQuery((PetscObject)dac0,"THI",(PetscObject*)&thi);
1272: if (!thi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot refine this DMDA, missing composed THI instance");
1273: if (nlevels > 1) {
1274: DMRefineHierarchy(dac0,nlevels-1,hierarchy);
1275: dac = hierarchy[nlevels-2];
1276: } else {
1277: dac = dac0;
1278: }
1279: DMDAGetInfo(dac,&dim, &N,&M,0, &n,&m,0, &dof,&s,0,0,0,&st);
1280: if (dim != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"This function can only refine 2D DMDAs");
1282: /* Creates a 3D DMDA with the same map-plane layout as the 2D one, with contiguous columns */
1283: 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);
1284: DMSetUp(daf);
1286: daf->ops->creatematrix = dac->ops->creatematrix;
1287: daf->ops->createinterpolation = dac->ops->createinterpolation;
1288: daf->ops->getcoloring = dac->ops->getcoloring;
1289: ddf = (DM_DA*)daf->data;
1290: ddc = (DM_DA*)dac->data;
1291: ddf->interptype = ddc->interptype;
1293: DMDASetFieldName(daf,0,"x-velocity");
1294: DMDASetFieldName(daf,1,"y-velocity");
1296: hierarchy[nlevels-1] = daf;
1297: return(0);
1298: }
1300: static PetscErrorCode DMCreateInterpolation_DA_THI(DM dac,DM daf,Mat *A,Vec *scale)
1301: {
1303: PetscInt dim;
1310: DMDAGetInfo(daf,&dim,0,0,0,0,0,0,0,0,0,0,0,0);
1311: if (dim == 2) {
1312: /* We are in the 2D problem and use normal DMDA interpolation */
1313: DMCreateInterpolation(dac,daf,A,scale);
1314: } else {
1315: PetscInt i,j,k,xs,ys,zs,xm,ym,zm,mx,my,mz,rstart,cstart;
1316: Mat B;
1318: DMDAGetInfo(daf,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
1319: DMDAGetCorners(daf,&zs,&ys,&xs,&zm,&ym,&xm);
1320: if (zs != 0) SETERRQ(PETSC_COMM_SELF,1,"unexpected");
1321: MatCreate(PetscObjectComm((PetscObject)daf),&B);
1322: MatSetSizes(B,xm*ym*zm,xm*ym,mx*my*mz,mx*my);
1324: MatSetType(B,MATAIJ);
1325: MatSeqAIJSetPreallocation(B,1,NULL);
1326: MatMPIAIJSetPreallocation(B,1,NULL,0,NULL);
1327: MatGetOwnershipRange(B,&rstart,NULL);
1328: MatGetOwnershipRangeColumn(B,&cstart,NULL);
1329: for (i=xs; i<xs+xm; i++) {
1330: for (j=ys; j<ys+ym; j++) {
1331: for (k=zs; k<zs+zm; k++) {
1332: PetscInt i2 = i*ym+j,i3 = i2*zm+k;
1333: PetscScalar val = ((k == 0 || k == mz-1) ? 0.5 : 1.) / (mz-1.); /* Integration using trapezoid rule */
1334: MatSetValue(B,cstart+i3,rstart+i2,val,INSERT_VALUES);
1335: }
1336: }
1337: }
1338: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1339: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1340: MatCreateMAIJ(B,sizeof(Node)/sizeof(PetscScalar),A);
1341: MatDestroy(&B);
1342: }
1343: return(0);
1344: }
1346: static PetscErrorCode DMCreateMatrix_THI_Tridiagonal(DM da,Mat *J)
1347: {
1348: PetscErrorCode ierr;
1349: Mat A;
1350: PetscInt xm,ym,zm,dim,dof = 2,starts[3],dims[3];
1351: ISLocalToGlobalMapping ltog;
1354: DMDAGetInfo(da,&dim, 0,0,0, 0,0,0, 0,0,0,0,0,0);
1355: if (dim != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected DMDA to be 3D");
1356: DMDAGetCorners(da,0,0,0,&zm,&ym,&xm);
1357: DMGetLocalToGlobalMapping(da,<og);
1358: MatCreate(PetscObjectComm((PetscObject)da),&A);
1359: MatSetSizes(A,dof*xm*ym*zm,dof*xm*ym*zm,PETSC_DETERMINE,PETSC_DETERMINE);
1360: MatSetType(A,da->mattype);
1361: MatSetFromOptions(A);
1362: MatSeqAIJSetPreallocation(A,3*2,NULL);
1363: MatMPIAIJSetPreallocation(A,3*2,NULL,0,NULL);
1364: MatSeqBAIJSetPreallocation(A,2,3,NULL);
1365: MatMPIBAIJSetPreallocation(A,2,3,NULL,0,NULL);
1366: MatSeqSBAIJSetPreallocation(A,2,2,NULL);
1367: MatMPISBAIJSetPreallocation(A,2,2,NULL,0,NULL);
1368: MatSetLocalToGlobalMapping(A,ltog,ltog);
1369: DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);
1370: MatSetStencil(A,dim,dims,starts,dof);
1371: *J = A;
1372: return(0);
1373: }
1375: static PetscErrorCode THIDAVecView_VTK_XML(THI thi,DM da,Vec X,const char filename[])
1376: {
1377: const PetscInt dof = 2;
1378: Units units = thi->units;
1379: MPI_Comm comm;
1380: PetscErrorCode ierr;
1381: PetscViewer viewer;
1382: PetscMPIInt rank,size,tag,nn,nmax;
1383: PetscInt mx,my,mz,r,range[6];
1384: const PetscScalar *x;
1387: PetscObjectGetComm((PetscObject)thi,&comm);
1388: DMDAGetInfo(da,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
1389: MPI_Comm_size(comm,&size);
1390: MPI_Comm_rank(comm,&rank);
1391: PetscViewerASCIIOpen(comm,filename,&viewer);
1392: PetscViewerASCIIPrintf(viewer,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");
1393: PetscViewerASCIIPrintf(viewer," <StructuredGrid WholeExtent=\"%d %D %d %D %d %D\">\n",0,mz-1,0,my-1,0,mx-1);
1395: DMDAGetCorners(da,range,range+1,range+2,range+3,range+4,range+5);
1396: PetscMPIIntCast(range[3]*range[4]*range[5]*dof,&nn);
1397: MPI_Reduce(&nn,&nmax,1,MPI_INT,MPI_MAX,0,comm);
1398: tag = ((PetscObject) viewer)->tag;
1399: VecGetArrayRead(X,&x);
1400: if (!rank) {
1401: PetscScalar *array;
1402: PetscMalloc1(nmax,&array);
1403: for (r=0; r<size; r++) {
1404: PetscInt i,j,k,xs,xm,ys,ym,zs,zm;
1405: const PetscScalar *ptr;
1406: MPI_Status status;
1407: if (r) {
1408: MPI_Recv(range,6,MPIU_INT,r,tag,comm,MPI_STATUS_IGNORE);
1409: }
1410: zs = range[0];ys = range[1];xs = range[2];zm = range[3];ym = range[4];xm = range[5];
1411: if (xm*ym*zm*dof > nmax) SETERRQ(PETSC_COMM_SELF,1,"should not happen");
1412: if (r) {
1413: MPI_Recv(array,nmax,MPIU_SCALAR,r,tag,comm,&status);
1414: MPI_Get_count(&status,MPIU_SCALAR,&nn);
1415: if (nn != xm*ym*zm*dof) SETERRQ(PETSC_COMM_SELF,1,"should not happen");
1416: ptr = array;
1417: } else ptr = x;
1418: PetscViewerASCIIPrintf(viewer," <Piece Extent=\"%D %D %D %D %D %D\">\n",zs,zs+zm-1,ys,ys+ym-1,xs,xs+xm-1);
1420: PetscViewerASCIIPrintf(viewer," <Points>\n");
1421: PetscViewerASCIIPrintf(viewer," <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n");
1422: for (i=xs; i<xs+xm; i++) {
1423: for (j=ys; j<ys+ym; j++) {
1424: for (k=zs; k<zs+zm; k++) {
1425: PrmNode p;
1426: PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my,zz;
1427: thi->initialize(thi,xx,yy,&p);
1428: zz = PetscRealPart(p.b) + PetscRealPart(p.h)*k/(mz-1);
1429: PetscViewerASCIIPrintf(viewer,"%f %f %f\n",(double)xx,(double)yy,(double)zz);
1430: }
1431: }
1432: }
1433: PetscViewerASCIIPrintf(viewer," </DataArray>\n");
1434: PetscViewerASCIIPrintf(viewer," </Points>\n");
1436: PetscViewerASCIIPrintf(viewer," <PointData>\n");
1437: PetscViewerASCIIPrintf(viewer," <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n");
1438: for (i=0; i<nn; i+=dof) {
1439: 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);
1440: }
1441: PetscViewerASCIIPrintf(viewer," </DataArray>\n");
1443: PetscViewerASCIIPrintf(viewer," <DataArray type=\"Int32\" Name=\"rank\" NumberOfComponents=\"1\" format=\"ascii\">\n");
1444: for (i=0; i<nn; i+=dof) {
1445: PetscViewerASCIIPrintf(viewer,"%D\n",r);
1446: }
1447: PetscViewerASCIIPrintf(viewer," </DataArray>\n");
1448: PetscViewerASCIIPrintf(viewer," </PointData>\n");
1450: PetscViewerASCIIPrintf(viewer," </Piece>\n");
1451: }
1452: PetscFree(array);
1453: } else {
1454: MPI_Send(range,6,MPIU_INT,0,tag,comm);
1455: MPI_Send((PetscScalar*)x,nn,MPIU_SCALAR,0,tag,comm);
1456: }
1457: VecRestoreArrayRead(X,&x);
1458: PetscViewerASCIIPrintf(viewer," </StructuredGrid>\n");
1459: PetscViewerASCIIPrintf(viewer,"</VTKFile>\n");
1460: PetscViewerDestroy(&viewer);
1461: return(0);
1462: }
1464: int main(int argc,char *argv[])
1465: {
1466: MPI_Comm comm;
1467: THI thi;
1469: DM da;
1470: SNES snes;
1472: PetscInitialize(&argc,&argv,0,help);
1473: comm = PETSC_COMM_WORLD;
1475: THICreate(comm,&thi);
1476: {
1477: PetscInt M = 3,N = 3,P = 2;
1478: PetscOptionsBegin(comm,NULL,"Grid resolution options","");
1479: {
1480: PetscOptionsInt("-M","Number of elements in x-direction on coarse level","",M,&M,NULL);
1481: N = M;
1482: PetscOptionsInt("-N","Number of elements in y-direction on coarse level (if different from M)","",N,&N,NULL);
1483: if (thi->coarse2d) {
1484: PetscOptionsInt("-zlevels","Number of elements in z-direction on fine level","",thi->zlevels,&thi->zlevels,NULL);
1485: } else {
1486: PetscOptionsInt("-P","Number of elements in z-direction on coarse level","",P,&P,NULL);
1487: }
1488: }
1489: PetscOptionsEnd();
1490: if (thi->coarse2d) {
1491: 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);
1492: DMSetFromOptions(da);
1493: DMSetUp(da);
1494: da->ops->refinehierarchy = DMRefineHierarchy_THI;
1495: da->ops->createinterpolation = DMCreateInterpolation_DA_THI;
1497: PetscObjectCompose((PetscObject)da,"THI",(PetscObject)thi);
1498: } else {
1499: 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);
1500: DMSetFromOptions(da);
1501: DMSetUp(da);
1502: }
1503: DMDASetFieldName(da,0,"x-velocity");
1504: DMDASetFieldName(da,1,"y-velocity");
1505: }
1506: THISetUpDM(thi,da);
1507: if (thi->tridiagonal) da->ops->creatematrix = DMCreateMatrix_THI_Tridiagonal;
1509: { /* Set the fine level matrix type if -da_refine */
1510: PetscInt rlevel,clevel;
1511: DMGetRefineLevel(da,&rlevel);
1512: DMGetCoarsenLevel(da,&clevel);
1513: if (rlevel - clevel > 0) {DMSetMatType(da,thi->mattype);}
1514: }
1516: DMDASNESSetFunctionLocal(da,ADD_VALUES,(DMDASNESFunction)THIFunctionLocal,thi);
1517: if (thi->tridiagonal) {
1518: DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)THIJacobianLocal_3D_Tridiagonal,thi);
1519: } else {
1520: DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)THIJacobianLocal_3D_Full,thi);
1521: }
1522: DMCoarsenHookAdd(da,DMCoarsenHook_THI,NULL,thi);
1523: DMRefineHookAdd(da,DMRefineHook_THI,NULL,thi);
1525: DMSetApplicationContext(da,thi);
1527: SNESCreate(comm,&snes);
1528: SNESSetDM(snes,da);
1529: DMDestroy(&da);
1530: SNESSetComputeInitialGuess(snes,THIInitial,NULL);
1531: SNESSetFromOptions(snes);
1533: SNESSolve(snes,NULL,NULL);
1535: THISolveStatistics(thi,snes,0,"Full");
1537: {
1538: PetscBool flg;
1539: char filename[PETSC_MAX_PATH_LEN] = "";
1540: PetscOptionsGetString(NULL,NULL,"-o",filename,sizeof(filename),&flg);
1541: if (flg) {
1542: Vec X;
1543: DM dm;
1544: SNESGetSolution(snes,&X);
1545: SNESGetDM(snes,&dm);
1546: THIDAVecView_VTK_XML(thi,dm,X,filename);
1547: }
1548: }
1550: DMDestroy(&da);
1551: SNESDestroy(&snes);
1552: THIDestroy(&thi);
1553: PetscFinalize();
1554: return ierr;
1555: }