xref: /petsc/src/snes/tutorials/ex28.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static const char help[] = "1D multiphysics prototype with analytic Jacobians to solve individual problems and a coupled problem.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /* Solve a PDE coupled to an algebraic system in 1D
4*c4762a1bSJed Brown  *
5*c4762a1bSJed Brown  * PDE (U):
6*c4762a1bSJed Brown  *     -(k u_x)_x = 1 on (0,1), subject to u(0) = 0, u(1) = 1
7*c4762a1bSJed Brown  * Algebraic (K):
8*c4762a1bSJed Brown  *     exp(k-1) + k = 1/(1/(1+u) + 1/(1+u_x^2))
9*c4762a1bSJed Brown  *
10*c4762a1bSJed Brown  * The discretization places k at staggered points, and a separate DMDA is used for each "physics".
11*c4762a1bSJed Brown  *
12*c4762a1bSJed Brown  * This example is a prototype for coupling in multi-physics problems, therefore residual evaluation and assembly for
13*c4762a1bSJed Brown  * each problem (referred to as U and K) are written separately.  This permits the same "physics" code to be used for
14*c4762a1bSJed Brown  * solving each uncoupled problem as well as the coupled system.  In particular, run with -problem_type 0 to solve only
15*c4762a1bSJed Brown  * problem U (with K fixed), -problem_type 1 to solve only K (with U fixed), and -problem_type 2 to solve both at once.
16*c4762a1bSJed Brown  *
17*c4762a1bSJed Brown  * In all cases, a fully-assembled analytic Jacobian is available, so the systems can be solved with a direct solve or
18*c4762a1bSJed Brown  * any other standard method.  Additionally, by running with
19*c4762a1bSJed Brown  *
20*c4762a1bSJed Brown  *   -pack_dm_mat_type nest
21*c4762a1bSJed Brown  *
22*c4762a1bSJed Brown  * The same code assembles a coupled matrix where each block is stored separately, which allows the use of PCFieldSplit
23*c4762a1bSJed Brown  * without copying values to extract submatrices.
24*c4762a1bSJed Brown  */
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown #include <petscsnes.h>
27*c4762a1bSJed Brown #include <petscdm.h>
28*c4762a1bSJed Brown #include <petscdmda.h>
29*c4762a1bSJed Brown #include <petscdmcomposite.h>
30*c4762a1bSJed Brown 
31*c4762a1bSJed Brown typedef struct _UserCtx *User;
32*c4762a1bSJed Brown struct _UserCtx {
33*c4762a1bSJed Brown   PetscInt ptype;
34*c4762a1bSJed Brown   DM       pack;
35*c4762a1bSJed Brown   Vec      Uloc,Kloc;
36*c4762a1bSJed Brown };
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal_U(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],PetscScalar f[])
39*c4762a1bSJed Brown {
40*c4762a1bSJed Brown   PetscReal hx = 1./info->mx;
41*c4762a1bSJed Brown   PetscInt  i;
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   PetscFunctionBeginUser;
44*c4762a1bSJed Brown   for (i=info->xs; i<info->xs+info->xm; i++) {
45*c4762a1bSJed Brown     if (i == 0) f[i] = 1./hx*u[i];
46*c4762a1bSJed Brown     else if (i == info->mx-1) f[i] = 1./hx*(u[i] - 1.0);
47*c4762a1bSJed Brown     else f[i] = hx*((k[i-1]*(u[i]-u[i-1]) - k[i]*(u[i+1]-u[i]))/(hx*hx) - 1.0);
48*c4762a1bSJed Brown   }
49*c4762a1bSJed Brown   PetscFunctionReturn(0);
50*c4762a1bSJed Brown }
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal_K(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],PetscScalar f[])
53*c4762a1bSJed Brown {
54*c4762a1bSJed Brown   PetscReal hx = 1./info->mx;
55*c4762a1bSJed Brown   PetscInt  i;
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown   PetscFunctionBeginUser;
58*c4762a1bSJed Brown   for (i=info->xs; i<info->xs+info->xm; i++) {
59*c4762a1bSJed Brown     const PetscScalar
60*c4762a1bSJed Brown       ubar  = 0.5*(u[i+1]+u[i]),
61*c4762a1bSJed Brown       gradu = (u[i+1]-u[i])/hx,
62*c4762a1bSJed Brown       g     = 1. + gradu*gradu,
63*c4762a1bSJed Brown       w     = 1./(1.+ubar) + 1./g;
64*c4762a1bSJed Brown     f[i] = hx*(PetscExpScalar(k[i]-1.0) + k[i] - 1./w);
65*c4762a1bSJed Brown   }
66*c4762a1bSJed Brown   PetscFunctionReturn(0);
67*c4762a1bSJed Brown }
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown static PetscErrorCode FormFunction_All(SNES snes,Vec X,Vec F,void *ctx)
70*c4762a1bSJed Brown {
71*c4762a1bSJed Brown   User           user = (User)ctx;
72*c4762a1bSJed Brown   DM             dau,dak;
73*c4762a1bSJed Brown   DMDALocalInfo  infou,infok;
74*c4762a1bSJed Brown   PetscScalar    *u,*k;
75*c4762a1bSJed Brown   PetscScalar    *fu,*fk;
76*c4762a1bSJed Brown   PetscErrorCode ierr;
77*c4762a1bSJed Brown   Vec            Uloc,Kloc,Fu,Fk;
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   PetscFunctionBeginUser;
80*c4762a1bSJed Brown   ierr = DMCompositeGetEntries(user->pack,&dau,&dak);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dau,&infou);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dak,&infok);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = DMCompositeGetLocalVectors(user->pack,&Uloc,&Kloc);CHKERRQ(ierr);
84*c4762a1bSJed Brown   switch (user->ptype) {
85*c4762a1bSJed Brown   case 0:
86*c4762a1bSJed Brown     ierr = DMGlobalToLocalBegin(dau,X,INSERT_VALUES,Uloc);CHKERRQ(ierr);
87*c4762a1bSJed Brown     ierr = DMGlobalToLocalEnd  (dau,X,INSERT_VALUES,Uloc);CHKERRQ(ierr);
88*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,Uloc,&u);CHKERRQ(ierr);
89*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,user->Kloc,&k);CHKERRQ(ierr);
90*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,F,&fu);CHKERRQ(ierr);
91*c4762a1bSJed Brown     ierr = FormFunctionLocal_U(user,&infou,u,k,fu);CHKERRQ(ierr);
92*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,F,&fu);CHKERRQ(ierr);
93*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,Uloc,&u);CHKERRQ(ierr);
94*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,user->Kloc,&k);CHKERRQ(ierr);
95*c4762a1bSJed Brown     break;
96*c4762a1bSJed Brown   case 1:
97*c4762a1bSJed Brown     ierr = DMGlobalToLocalBegin(dak,X,INSERT_VALUES,Kloc);CHKERRQ(ierr);
98*c4762a1bSJed Brown     ierr = DMGlobalToLocalEnd  (dak,X,INSERT_VALUES,Kloc);CHKERRQ(ierr);
99*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,user->Uloc,&u);CHKERRQ(ierr);
100*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,Kloc,&k);CHKERRQ(ierr);
101*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,F,&fk);CHKERRQ(ierr);
102*c4762a1bSJed Brown     ierr = FormFunctionLocal_K(user,&infok,u,k,fk);CHKERRQ(ierr);
103*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,F,&fk);CHKERRQ(ierr);
104*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,user->Uloc,&u);CHKERRQ(ierr);
105*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,Kloc,&k);CHKERRQ(ierr);
106*c4762a1bSJed Brown     break;
107*c4762a1bSJed Brown   case 2:
108*c4762a1bSJed Brown     ierr = DMCompositeScatter(user->pack,X,Uloc,Kloc);CHKERRQ(ierr);
109*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,Uloc,&u);CHKERRQ(ierr);
110*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,Kloc,&k);CHKERRQ(ierr);
111*c4762a1bSJed Brown     ierr = DMCompositeGetAccess(user->pack,F,&Fu,&Fk);CHKERRQ(ierr);
112*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,Fu,&fu);CHKERRQ(ierr);
113*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,Fk,&fk);CHKERRQ(ierr);
114*c4762a1bSJed Brown     ierr = FormFunctionLocal_U(user,&infou,u,k,fu);CHKERRQ(ierr);
115*c4762a1bSJed Brown     ierr = FormFunctionLocal_K(user,&infok,u,k,fk);CHKERRQ(ierr);
116*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,Fu,&fu);CHKERRQ(ierr);
117*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,Fk,&fk);CHKERRQ(ierr);
118*c4762a1bSJed Brown     ierr = DMCompositeRestoreAccess(user->pack,F,&Fu,&Fk);CHKERRQ(ierr);
119*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,Uloc,&u);CHKERRQ(ierr);
120*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,Kloc,&k);CHKERRQ(ierr);
121*c4762a1bSJed Brown     break;
122*c4762a1bSJed Brown   }
123*c4762a1bSJed Brown   ierr = DMCompositeRestoreLocalVectors(user->pack,&Uloc,&Kloc);CHKERRQ(ierr);
124*c4762a1bSJed Brown   PetscFunctionReturn(0);
125*c4762a1bSJed Brown }
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal_U(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],Mat Buu)
128*c4762a1bSJed Brown {
129*c4762a1bSJed Brown   PetscReal      hx = 1./info->mx;
130*c4762a1bSJed Brown   PetscErrorCode ierr;
131*c4762a1bSJed Brown   PetscInt       i;
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   PetscFunctionBeginUser;
134*c4762a1bSJed Brown   for (i=info->xs; i<info->xs+info->xm; i++) {
135*c4762a1bSJed Brown     PetscInt    row = i-info->gxs,cols[] = {row-1,row,row+1};
136*c4762a1bSJed Brown     PetscScalar val = 1./hx;
137*c4762a1bSJed Brown     if (i == 0) {
138*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Buu,1,&row,1,&row,&val,INSERT_VALUES);CHKERRQ(ierr);
139*c4762a1bSJed Brown     } else if (i == info->mx-1) {
140*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Buu,1,&row,1,&row,&val,INSERT_VALUES);CHKERRQ(ierr);
141*c4762a1bSJed Brown     } else {
142*c4762a1bSJed Brown       PetscScalar vals[] = {-k[i-1]/hx,(k[i-1]+k[i])/hx,-k[i]/hx};
143*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Buu,1,&row,3,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
144*c4762a1bSJed Brown     }
145*c4762a1bSJed Brown   }
146*c4762a1bSJed Brown   PetscFunctionReturn(0);
147*c4762a1bSJed Brown }
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal_K(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],Mat Bkk)
150*c4762a1bSJed Brown {
151*c4762a1bSJed Brown   PetscReal      hx = 1./info->mx;
152*c4762a1bSJed Brown   PetscErrorCode ierr;
153*c4762a1bSJed Brown   PetscInt       i;
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown   PetscFunctionBeginUser;
156*c4762a1bSJed Brown   for (i=info->xs; i<info->xs+info->xm; i++) {
157*c4762a1bSJed Brown     PetscInt    row    = i-info->gxs;
158*c4762a1bSJed Brown     PetscScalar vals[] = {hx*(PetscExpScalar(k[i]-1.)+ (PetscScalar)1.)};
159*c4762a1bSJed Brown     ierr = MatSetValuesLocal(Bkk,1,&row,1,&row,vals,INSERT_VALUES);CHKERRQ(ierr);
160*c4762a1bSJed Brown   }
161*c4762a1bSJed Brown   PetscFunctionReturn(0);
162*c4762a1bSJed Brown }
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal_UK(User user,DMDALocalInfo *info,DMDALocalInfo *infok,const PetscScalar u[],const PetscScalar k[],Mat Buk)
165*c4762a1bSJed Brown {
166*c4762a1bSJed Brown   PetscReal      hx = 1./info->mx;
167*c4762a1bSJed Brown   PetscErrorCode ierr;
168*c4762a1bSJed Brown   PetscInt       i;
169*c4762a1bSJed Brown   PetscInt       row,cols[2];
170*c4762a1bSJed Brown   PetscScalar    vals[2];
171*c4762a1bSJed Brown 
172*c4762a1bSJed Brown   PetscFunctionBeginUser;
173*c4762a1bSJed Brown   if (!Buk) PetscFunctionReturn(0); /* Not assembling this block */
174*c4762a1bSJed Brown   for (i=info->xs; i<info->xs+info->xm; i++) {
175*c4762a1bSJed Brown     if (i == 0 || i == info->mx-1) continue;
176*c4762a1bSJed Brown     row     = i-info->gxs;
177*c4762a1bSJed Brown     cols[0] = i-1-infok->gxs;  vals[0] = (u[i]-u[i-1])/hx;
178*c4762a1bSJed Brown     cols[1] = i-infok->gxs;    vals[1] = (u[i]-u[i+1])/hx;
179*c4762a1bSJed Brown     ierr    = MatSetValuesLocal(Buk,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
180*c4762a1bSJed Brown   }
181*c4762a1bSJed Brown   PetscFunctionReturn(0);
182*c4762a1bSJed Brown }
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal_KU(User user,DMDALocalInfo *info,DMDALocalInfo *infok,const PetscScalar u[],const PetscScalar k[],Mat Bku)
185*c4762a1bSJed Brown {
186*c4762a1bSJed Brown   PetscErrorCode ierr;
187*c4762a1bSJed Brown   PetscInt       i;
188*c4762a1bSJed Brown   PetscReal      hx = 1./(info->mx-1);
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown   PetscFunctionBeginUser;
191*c4762a1bSJed Brown   if (!Bku) PetscFunctionReturn(0); /* Not assembling this block */
192*c4762a1bSJed Brown   for (i=infok->xs; i<infok->xs+infok->xm; i++) {
193*c4762a1bSJed Brown     PetscInt    row = i-infok->gxs,cols[2];
194*c4762a1bSJed Brown     PetscScalar vals[2];
195*c4762a1bSJed Brown     const PetscScalar
196*c4762a1bSJed Brown       ubar     = 0.5*(u[i]+u[i+1]),
197*c4762a1bSJed Brown       ubar_L   = 0.5,
198*c4762a1bSJed Brown       ubar_R   = 0.5,
199*c4762a1bSJed Brown       gradu    = (u[i+1]-u[i])/hx,
200*c4762a1bSJed Brown       gradu_L  = -1./hx,
201*c4762a1bSJed Brown       gradu_R  = 1./hx,
202*c4762a1bSJed Brown       g        = 1. + PetscSqr(gradu),
203*c4762a1bSJed Brown       g_gradu  = 2.*gradu,
204*c4762a1bSJed Brown       w        = 1./(1.+ubar) + 1./g,
205*c4762a1bSJed Brown       w_ubar   = -1./PetscSqr(1.+ubar),
206*c4762a1bSJed Brown       w_gradu  = -g_gradu/PetscSqr(g),
207*c4762a1bSJed Brown       iw       = 1./w,
208*c4762a1bSJed Brown       iw_ubar  = -w_ubar * PetscSqr(iw),
209*c4762a1bSJed Brown       iw_gradu = -w_gradu * PetscSqr(iw);
210*c4762a1bSJed Brown     cols[0] = i-info->gxs;         vals[0] = -hx*(iw_ubar*ubar_L + iw_gradu*gradu_L);
211*c4762a1bSJed Brown     cols[1] = i+1-info->gxs;       vals[1] = -hx*(iw_ubar*ubar_R + iw_gradu*gradu_R);
212*c4762a1bSJed Brown     ierr    = MatSetValuesLocal(Bku,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
213*c4762a1bSJed Brown   }
214*c4762a1bSJed Brown   PetscFunctionReturn(0);
215*c4762a1bSJed Brown }
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown static PetscErrorCode FormJacobian_All(SNES snes,Vec X,Mat J,Mat B,void *ctx)
218*c4762a1bSJed Brown {
219*c4762a1bSJed Brown   User           user = (User)ctx;
220*c4762a1bSJed Brown   DM             dau,dak;
221*c4762a1bSJed Brown   DMDALocalInfo  infou,infok;
222*c4762a1bSJed Brown   PetscScalar    *u,*k;
223*c4762a1bSJed Brown   PetscErrorCode ierr;
224*c4762a1bSJed Brown   Vec            Uloc,Kloc;
225*c4762a1bSJed Brown 
226*c4762a1bSJed Brown   PetscFunctionBeginUser;
227*c4762a1bSJed Brown   ierr = DMCompositeGetEntries(user->pack,&dau,&dak);CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dau,&infou);CHKERRQ(ierr);
229*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dak,&infok);CHKERRQ(ierr);
230*c4762a1bSJed Brown   ierr = DMCompositeGetLocalVectors(user->pack,&Uloc,&Kloc);CHKERRQ(ierr);
231*c4762a1bSJed Brown   switch (user->ptype) {
232*c4762a1bSJed Brown   case 0:
233*c4762a1bSJed Brown     ierr = DMGlobalToLocalBegin(dau,X,INSERT_VALUES,Uloc);CHKERRQ(ierr);
234*c4762a1bSJed Brown     ierr = DMGlobalToLocalEnd  (dau,X,INSERT_VALUES,Uloc);CHKERRQ(ierr);
235*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,Uloc,&u);CHKERRQ(ierr);
236*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,user->Kloc,&k);CHKERRQ(ierr);
237*c4762a1bSJed Brown     ierr = FormJacobianLocal_U(user,&infou,u,k,B);CHKERRQ(ierr);
238*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,Uloc,&u);CHKERRQ(ierr);
239*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,user->Kloc,&k);CHKERRQ(ierr);
240*c4762a1bSJed Brown     break;
241*c4762a1bSJed Brown   case 1:
242*c4762a1bSJed Brown     ierr = DMGlobalToLocalBegin(dak,X,INSERT_VALUES,Kloc);CHKERRQ(ierr);
243*c4762a1bSJed Brown     ierr = DMGlobalToLocalEnd  (dak,X,INSERT_VALUES,Kloc);CHKERRQ(ierr);
244*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,user->Uloc,&u);CHKERRQ(ierr);
245*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,Kloc,&k);CHKERRQ(ierr);
246*c4762a1bSJed Brown     ierr = FormJacobianLocal_K(user,&infok,u,k,B);CHKERRQ(ierr);
247*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,user->Uloc,&u);CHKERRQ(ierr);
248*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,Kloc,&k);CHKERRQ(ierr);
249*c4762a1bSJed Brown     break;
250*c4762a1bSJed Brown   case 2: {
251*c4762a1bSJed Brown     Mat       Buu,Buk,Bku,Bkk;
252*c4762a1bSJed Brown     PetscBool nest;
253*c4762a1bSJed Brown     IS        *is;
254*c4762a1bSJed Brown     ierr = DMCompositeScatter(user->pack,X,Uloc,Kloc);CHKERRQ(ierr);
255*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dau,Uloc,&u);CHKERRQ(ierr);
256*c4762a1bSJed Brown     ierr = DMDAVecGetArray(dak,Kloc,&k);CHKERRQ(ierr);
257*c4762a1bSJed Brown     ierr = DMCompositeGetLocalISs(user->pack,&is);CHKERRQ(ierr);
258*c4762a1bSJed Brown     ierr = MatGetLocalSubMatrix(B,is[0],is[0],&Buu);CHKERRQ(ierr);
259*c4762a1bSJed Brown     ierr = MatGetLocalSubMatrix(B,is[0],is[1],&Buk);CHKERRQ(ierr);
260*c4762a1bSJed Brown     ierr = MatGetLocalSubMatrix(B,is[1],is[0],&Bku);CHKERRQ(ierr);
261*c4762a1bSJed Brown     ierr = MatGetLocalSubMatrix(B,is[1],is[1],&Bkk);CHKERRQ(ierr);
262*c4762a1bSJed Brown     ierr = FormJacobianLocal_U(user,&infou,u,k,Buu);CHKERRQ(ierr);
263*c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)B,MATNEST,&nest);CHKERRQ(ierr);
264*c4762a1bSJed Brown     if (!nest) {
265*c4762a1bSJed Brown       /*
266*c4762a1bSJed Brown          DMCreateMatrix_Composite()  for a nested matrix does not generate off-block matrices that one can call MatSetValuesLocal() on, it just creates dummy
267*c4762a1bSJed Brown          matrices with no entries; there cannot insert values into them. Commit b6480e041dd2293a65f96222772d68cdb4ed6306
268*c4762a1bSJed Brown          changed Mat_Nest() from returning NULL pointers for these submatrices to dummy matrices because PCFIELDSPLIT could not
269*c4762a1bSJed Brown          handle the returned null matrices.
270*c4762a1bSJed Brown       */
271*c4762a1bSJed Brown       ierr = FormJacobianLocal_UK(user,&infou,&infok,u,k,Buk);CHKERRQ(ierr);
272*c4762a1bSJed Brown       ierr = FormJacobianLocal_KU(user,&infou,&infok,u,k,Bku);CHKERRQ(ierr);
273*c4762a1bSJed Brown     }
274*c4762a1bSJed Brown     ierr = FormJacobianLocal_K(user,&infok,u,k,Bkk);CHKERRQ(ierr);
275*c4762a1bSJed Brown     ierr = MatRestoreLocalSubMatrix(B,is[0],is[0],&Buu);CHKERRQ(ierr);
276*c4762a1bSJed Brown     ierr = MatRestoreLocalSubMatrix(B,is[0],is[1],&Buk);CHKERRQ(ierr);
277*c4762a1bSJed Brown     ierr = MatRestoreLocalSubMatrix(B,is[1],is[0],&Bku);CHKERRQ(ierr);
278*c4762a1bSJed Brown     ierr = MatRestoreLocalSubMatrix(B,is[1],is[1],&Bkk);CHKERRQ(ierr);
279*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dau,Uloc,&u);CHKERRQ(ierr);
280*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dak,Kloc,&k);CHKERRQ(ierr);
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown     ierr = ISDestroy(&is[0]);CHKERRQ(ierr);
283*c4762a1bSJed Brown     ierr = ISDestroy(&is[1]);CHKERRQ(ierr);
284*c4762a1bSJed Brown     ierr = PetscFree(is);CHKERRQ(ierr);
285*c4762a1bSJed Brown   } break;
286*c4762a1bSJed Brown   }
287*c4762a1bSJed Brown   ierr = DMCompositeRestoreLocalVectors(user->pack,&Uloc,&Kloc);CHKERRQ(ierr);
288*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
289*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
290*c4762a1bSJed Brown   if (J != B) {
291*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
292*c4762a1bSJed Brown     ierr = MatAssemblyEnd  (J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293*c4762a1bSJed Brown   }
294*c4762a1bSJed Brown   PetscFunctionReturn(0);
295*c4762a1bSJed Brown }
296*c4762a1bSJed Brown 
297*c4762a1bSJed Brown static PetscErrorCode FormInitial_Coupled(User user,Vec X)
298*c4762a1bSJed Brown {
299*c4762a1bSJed Brown   PetscErrorCode ierr;
300*c4762a1bSJed Brown   DM             dau,dak;
301*c4762a1bSJed Brown   DMDALocalInfo  infou,infok;
302*c4762a1bSJed Brown   Vec            Xu,Xk;
303*c4762a1bSJed Brown   PetscScalar    *u,*k,hx;
304*c4762a1bSJed Brown   PetscInt       i;
305*c4762a1bSJed Brown 
306*c4762a1bSJed Brown   PetscFunctionBeginUser;
307*c4762a1bSJed Brown   ierr = DMCompositeGetEntries(user->pack,&dau,&dak);CHKERRQ(ierr);
308*c4762a1bSJed Brown   ierr = DMCompositeGetAccess(user->pack,X,&Xu,&Xk);CHKERRQ(ierr);
309*c4762a1bSJed Brown   ierr = DMDAVecGetArray(dau,Xu,&u);CHKERRQ(ierr);
310*c4762a1bSJed Brown   ierr = DMDAVecGetArray(dak,Xk,&k);CHKERRQ(ierr);
311*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dau,&infou);CHKERRQ(ierr);
312*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dak,&infok);CHKERRQ(ierr);
313*c4762a1bSJed Brown   hx   = 1./(infok.mx);
314*c4762a1bSJed Brown   for (i=infou.xs; i<infou.xs+infou.xm; i++) u[i] = (PetscScalar)i*hx * (1.-(PetscScalar)i*hx);
315*c4762a1bSJed Brown   for (i=infok.xs; i<infok.xs+infok.xm; i++) k[i] = 1.0 + 0.5*PetscSinScalar(2*PETSC_PI*i*hx);
316*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(dau,Xu,&u);CHKERRQ(ierr);
317*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(dak,Xk,&k);CHKERRQ(ierr);
318*c4762a1bSJed Brown   ierr = DMCompositeRestoreAccess(user->pack,X,&Xu,&Xk);CHKERRQ(ierr);
319*c4762a1bSJed Brown   ierr = DMCompositeScatter(user->pack,X,user->Uloc,user->Kloc);CHKERRQ(ierr);
320*c4762a1bSJed Brown   PetscFunctionReturn(0);
321*c4762a1bSJed Brown }
322*c4762a1bSJed Brown 
323*c4762a1bSJed Brown int main(int argc, char *argv[])
324*c4762a1bSJed Brown {
325*c4762a1bSJed Brown   PetscErrorCode ierr;
326*c4762a1bSJed Brown   DM             dau,dak,pack;
327*c4762a1bSJed Brown   const PetscInt *lxu;
328*c4762a1bSJed Brown   PetscInt       *lxk,m,sizes;
329*c4762a1bSJed Brown   User           user;
330*c4762a1bSJed Brown   SNES           snes;
331*c4762a1bSJed Brown   Vec            X,F,Xu,Xk,Fu,Fk;
332*c4762a1bSJed Brown   Mat            B;
333*c4762a1bSJed Brown   IS             *isg;
334*c4762a1bSJed Brown   PetscBool      pass_dm;
335*c4762a1bSJed Brown 
336*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
337*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,10,1,1,NULL,&dau);CHKERRQ(ierr);
338*c4762a1bSJed Brown   ierr = DMSetOptionsPrefix(dau,"u_");CHKERRQ(ierr);
339*c4762a1bSJed Brown   ierr = DMSetFromOptions(dau);CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = DMSetUp(dau);CHKERRQ(ierr);
341*c4762a1bSJed Brown   ierr = DMDAGetOwnershipRanges(dau,&lxu,0,0);CHKERRQ(ierr);
342*c4762a1bSJed Brown   ierr = DMDAGetInfo(dau,0, &m,0,0, &sizes,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
343*c4762a1bSJed Brown   ierr = PetscMalloc1(sizes,&lxk);CHKERRQ(ierr);
344*c4762a1bSJed Brown   ierr = PetscArraycpy(lxk,lxu,sizes);CHKERRQ(ierr);
345*c4762a1bSJed Brown   lxk[0]--;
346*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,m-1,1,1,lxk,&dak);CHKERRQ(ierr);
347*c4762a1bSJed Brown   ierr = DMSetOptionsPrefix(dak,"k_");CHKERRQ(ierr);
348*c4762a1bSJed Brown   ierr = DMSetFromOptions(dak);CHKERRQ(ierr);
349*c4762a1bSJed Brown   ierr = DMSetUp(dak);CHKERRQ(ierr);
350*c4762a1bSJed Brown   ierr = PetscFree(lxk);CHKERRQ(ierr);
351*c4762a1bSJed Brown 
352*c4762a1bSJed Brown   ierr = DMCompositeCreate(PETSC_COMM_WORLD,&pack);CHKERRQ(ierr);
353*c4762a1bSJed Brown   ierr = DMSetOptionsPrefix(pack,"pack_");CHKERRQ(ierr);
354*c4762a1bSJed Brown   ierr = DMCompositeAddDM(pack,dau);CHKERRQ(ierr);
355*c4762a1bSJed Brown   ierr = DMCompositeAddDM(pack,dak);CHKERRQ(ierr);
356*c4762a1bSJed Brown   ierr = DMDASetFieldName(dau,0,"u");CHKERRQ(ierr);
357*c4762a1bSJed Brown   ierr = DMDASetFieldName(dak,0,"k");CHKERRQ(ierr);
358*c4762a1bSJed Brown   ierr = DMSetFromOptions(pack);CHKERRQ(ierr);
359*c4762a1bSJed Brown 
360*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(pack,&X);CHKERRQ(ierr);
361*c4762a1bSJed Brown   ierr = VecDuplicate(X,&F);CHKERRQ(ierr);
362*c4762a1bSJed Brown 
363*c4762a1bSJed Brown   ierr = PetscNew(&user);CHKERRQ(ierr);
364*c4762a1bSJed Brown 
365*c4762a1bSJed Brown   user->pack = pack;
366*c4762a1bSJed Brown 
367*c4762a1bSJed Brown   ierr = DMCompositeGetGlobalISs(pack,&isg);CHKERRQ(ierr);
368*c4762a1bSJed Brown   ierr = DMCompositeGetLocalVectors(pack,&user->Uloc,&user->Kloc);CHKERRQ(ierr);
369*c4762a1bSJed Brown   ierr = DMCompositeScatter(pack,X,user->Uloc,user->Kloc);CHKERRQ(ierr);
370*c4762a1bSJed Brown 
371*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Coupled problem options","SNES");CHKERRQ(ierr);
372*c4762a1bSJed Brown   {
373*c4762a1bSJed Brown     user->ptype = 0; pass_dm = PETSC_TRUE;
374*c4762a1bSJed Brown 
375*c4762a1bSJed Brown     ierr = PetscOptionsInt("-problem_type","0: solve for u only, 1: solve for k only, 2: solve for both",0,user->ptype,&user->ptype,NULL);CHKERRQ(ierr);
376*c4762a1bSJed Brown     ierr = PetscOptionsBool("-pass_dm","Pass the packed DM to SNES to use when determining splits and forward into splits",0,pass_dm,&pass_dm,NULL);CHKERRQ(ierr);
377*c4762a1bSJed Brown   }
378*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
379*c4762a1bSJed Brown 
380*c4762a1bSJed Brown   ierr = FormInitial_Coupled(user,X);CHKERRQ(ierr);
381*c4762a1bSJed Brown 
382*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
383*c4762a1bSJed Brown   switch (user->ptype) {
384*c4762a1bSJed Brown   case 0:
385*c4762a1bSJed Brown     ierr = DMCompositeGetAccess(pack,X,&Xu,0);CHKERRQ(ierr);
386*c4762a1bSJed Brown     ierr = DMCompositeGetAccess(pack,F,&Fu,0);CHKERRQ(ierr);
387*c4762a1bSJed Brown     ierr = DMCreateMatrix(dau,&B);CHKERRQ(ierr);
388*c4762a1bSJed Brown     ierr = SNESSetFunction(snes,Fu,FormFunction_All,user);CHKERRQ(ierr);
389*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,B,B,FormJacobian_All,user);CHKERRQ(ierr);
390*c4762a1bSJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
391*c4762a1bSJed Brown     ierr = SNESSetDM(snes,dau);CHKERRQ(ierr);
392*c4762a1bSJed Brown     ierr = SNESSolve(snes,NULL,Xu);CHKERRQ(ierr);
393*c4762a1bSJed Brown     ierr = DMCompositeRestoreAccess(pack,X,&Xu,0);CHKERRQ(ierr);
394*c4762a1bSJed Brown     ierr = DMCompositeRestoreAccess(pack,F,&Fu,0);CHKERRQ(ierr);
395*c4762a1bSJed Brown     break;
396*c4762a1bSJed Brown   case 1:
397*c4762a1bSJed Brown     ierr = DMCompositeGetAccess(pack,X,0,&Xk);CHKERRQ(ierr);
398*c4762a1bSJed Brown     ierr = DMCompositeGetAccess(pack,F,0,&Fk);CHKERRQ(ierr);
399*c4762a1bSJed Brown     ierr = DMCreateMatrix(dak,&B);CHKERRQ(ierr);
400*c4762a1bSJed Brown     ierr = SNESSetFunction(snes,Fk,FormFunction_All,user);CHKERRQ(ierr);
401*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,B,B,FormJacobian_All,user);CHKERRQ(ierr);
402*c4762a1bSJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
403*c4762a1bSJed Brown     ierr = SNESSetDM(snes,dak);CHKERRQ(ierr);
404*c4762a1bSJed Brown     ierr = SNESSolve(snes,NULL,Xk);CHKERRQ(ierr);
405*c4762a1bSJed Brown     ierr = DMCompositeRestoreAccess(pack,X,0,&Xk);CHKERRQ(ierr);
406*c4762a1bSJed Brown     ierr = DMCompositeRestoreAccess(pack,F,0,&Fk);CHKERRQ(ierr);
407*c4762a1bSJed Brown     break;
408*c4762a1bSJed Brown   case 2:
409*c4762a1bSJed Brown     ierr = DMCreateMatrix(pack,&B);CHKERRQ(ierr);
410*c4762a1bSJed Brown     /* This example does not correctly allocate off-diagonal blocks. These options allows new nonzeros (slow). */
411*c4762a1bSJed Brown     ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
412*c4762a1bSJed Brown     ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
413*c4762a1bSJed Brown     ierr = SNESSetFunction(snes,F,FormFunction_All,user);CHKERRQ(ierr);
414*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,B,B,FormJacobian_All,user);CHKERRQ(ierr);
415*c4762a1bSJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
416*c4762a1bSJed Brown     if (!pass_dm) {             /* Manually provide index sets and names for the splits */
417*c4762a1bSJed Brown       KSP ksp;
418*c4762a1bSJed Brown       PC  pc;
419*c4762a1bSJed Brown       ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
420*c4762a1bSJed Brown       ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
421*c4762a1bSJed Brown       ierr = PCFieldSplitSetIS(pc,"u",isg[0]);CHKERRQ(ierr);
422*c4762a1bSJed Brown       ierr = PCFieldSplitSetIS(pc,"k",isg[1]);CHKERRQ(ierr);
423*c4762a1bSJed Brown     } else {
424*c4762a1bSJed Brown       /* The same names come from the options prefix for dau and dak. This option can support geometric multigrid inside
425*c4762a1bSJed Brown        * of splits, but it requires using a DM (perhaps your own implementation). */
426*c4762a1bSJed Brown       ierr = SNESSetDM(snes,pack);CHKERRQ(ierr);
427*c4762a1bSJed Brown     }
428*c4762a1bSJed Brown     ierr = SNESSolve(snes,NULL,X);CHKERRQ(ierr);
429*c4762a1bSJed Brown     break;
430*c4762a1bSJed Brown   }
431*c4762a1bSJed Brown   ierr = VecViewFromOptions(X,NULL,"-view_sol");CHKERRQ(ierr);
432*c4762a1bSJed Brown 
433*c4762a1bSJed Brown   if (0) {
434*c4762a1bSJed Brown     PetscInt  col      = 0;
435*c4762a1bSJed Brown     PetscBool mult_dup = PETSC_FALSE,view_dup = PETSC_FALSE;
436*c4762a1bSJed Brown     Mat       D;
437*c4762a1bSJed Brown     Vec       Y;
438*c4762a1bSJed Brown 
439*c4762a1bSJed Brown     ierr = PetscOptionsGetInt(NULL,0,"-col",&col,0);CHKERRQ(ierr);
440*c4762a1bSJed Brown     ierr = PetscOptionsGetBool(NULL,0,"-mult_dup",&mult_dup,0);CHKERRQ(ierr);
441*c4762a1bSJed Brown     ierr = PetscOptionsGetBool(NULL,0,"-view_dup",&view_dup,0);CHKERRQ(ierr);
442*c4762a1bSJed Brown 
443*c4762a1bSJed Brown     ierr = VecDuplicate(X,&Y);CHKERRQ(ierr);
444*c4762a1bSJed Brown     /* ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
445*c4762a1bSJed Brown     /* ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */
446*c4762a1bSJed Brown     ierr = MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr);
447*c4762a1bSJed Brown     ierr = VecZeroEntries(X);CHKERRQ(ierr);
448*c4762a1bSJed Brown     ierr = VecSetValue(X,col,1.0,INSERT_VALUES);CHKERRQ(ierr);
449*c4762a1bSJed Brown     ierr = VecAssemblyBegin(X);CHKERRQ(ierr);
450*c4762a1bSJed Brown     ierr = VecAssemblyEnd(X);CHKERRQ(ierr);
451*c4762a1bSJed Brown     ierr = MatMult(mult_dup ? D : B,X,Y);CHKERRQ(ierr);
452*c4762a1bSJed Brown     ierr = MatView(view_dup ? D : B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
453*c4762a1bSJed Brown     /* ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
454*c4762a1bSJed Brown     ierr = VecView(Y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
455*c4762a1bSJed Brown     ierr = MatDestroy(&D);CHKERRQ(ierr);
456*c4762a1bSJed Brown     ierr = VecDestroy(&Y);CHKERRQ(ierr);
457*c4762a1bSJed Brown   }
458*c4762a1bSJed Brown 
459*c4762a1bSJed Brown   ierr = DMCompositeRestoreLocalVectors(pack,&user->Uloc,&user->Kloc);CHKERRQ(ierr);
460*c4762a1bSJed Brown   ierr = PetscFree(user);CHKERRQ(ierr);
461*c4762a1bSJed Brown 
462*c4762a1bSJed Brown   ierr = ISDestroy(&isg[0]);CHKERRQ(ierr);
463*c4762a1bSJed Brown   ierr = ISDestroy(&isg[1]);CHKERRQ(ierr);
464*c4762a1bSJed Brown   ierr = PetscFree(isg);CHKERRQ(ierr);
465*c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
466*c4762a1bSJed Brown   ierr = VecDestroy(&F);CHKERRQ(ierr);
467*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
468*c4762a1bSJed Brown   ierr = DMDestroy(&dau);CHKERRQ(ierr);
469*c4762a1bSJed Brown   ierr = DMDestroy(&dak);CHKERRQ(ierr);
470*c4762a1bSJed Brown   ierr = DMDestroy(&pack);CHKERRQ(ierr);
471*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
472*c4762a1bSJed Brown   ierr = PetscFinalize();
473*c4762a1bSJed Brown   return ierr;
474*c4762a1bSJed Brown }
475*c4762a1bSJed Brown 
476*c4762a1bSJed Brown /*TEST
477*c4762a1bSJed Brown 
478*c4762a1bSJed Brown    build:
479*c4762a1bSJed Brown       requires: c99
480*c4762a1bSJed Brown 
481*c4762a1bSJed Brown    test:
482*c4762a1bSJed Brown       suffix: 0
483*c4762a1bSJed Brown       nsize: 3
484*c4762a1bSJed Brown       args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 0
485*c4762a1bSJed Brown 
486*c4762a1bSJed Brown    test:
487*c4762a1bSJed Brown       suffix: 1
488*c4762a1bSJed Brown       nsize: 3
489*c4762a1bSJed Brown       args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 1
490*c4762a1bSJed Brown 
491*c4762a1bSJed Brown    test:
492*c4762a1bSJed Brown       suffix: 2
493*c4762a1bSJed Brown       nsize: 3
494*c4762a1bSJed Brown       args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 2
495*c4762a1bSJed Brown 
496*c4762a1bSJed Brown    test:
497*c4762a1bSJed Brown       suffix: 3
498*c4762a1bSJed Brown       nsize: 3
499*c4762a1bSJed Brown       args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -ksp_monitor_short -problem_type 2 -snes_mf_operator -pack_dm_mat_type {{aij nest}} -pc_type fieldsplit -pc_fieldsplit_dm_splits -pc_fieldsplit_type additive -fieldsplit_u_ksp_type gmres -fieldsplit_k_pc_type jacobi
500*c4762a1bSJed Brown 
501*c4762a1bSJed Brown    test:
502*c4762a1bSJed Brown       suffix: 4
503*c4762a1bSJed Brown       nsize: 6
504*c4762a1bSJed Brown       args: -u_da_grid_x 257 -snes_converged_reason -snes_monitor_short -ksp_monitor_short -problem_type 2 -snes_mf_operator -pack_dm_mat_type aij -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_u_ksp_type gmres -fieldsplit_u_ksp_pc_side right -fieldsplit_u_pc_type mg -fieldsplit_u_pc_mg_levels 4 -fieldsplit_u_mg_levels_ksp_type richardson -fieldsplit_u_mg_levels_ksp_max_it 1 -fieldsplit_u_mg_levels_pc_type sor -fieldsplit_u_pc_mg_galerkin pmat -fieldsplit_u_ksp_converged_reason -fieldsplit_k_pc_type jacobi
505*c4762a1bSJed Brown       requires: !single
506*c4762a1bSJed Brown 
507*c4762a1bSJed Brown    test:
508*c4762a1bSJed Brown       suffix: glvis_composite_da_1d
509*c4762a1bSJed Brown       args: -u_da_grid_x 20 -problem_type 0 -snes_monitor_solution glvis:
510*c4762a1bSJed Brown 
511*c4762a1bSJed Brown TEST*/
512