xref: /libCEED/examples/solids/src/matops.c (revision 5754ecac3b7d1ff97b39b25dc78c06350f2c900d)
1ccaff030SJeremy L Thompson /// @file
2ccaff030SJeremy L Thompson /// Matrix shell operations for solid mechanics example using PETSc
3ccaff030SJeremy L Thompson 
4*5754ecacSJeremy L Thompson #include "../include/matops.h"
5*5754ecacSJeremy L Thompson #include "../include/structs.h"
6*5754ecacSJeremy L Thompson #include "../include/utils.h"
7ccaff030SJeremy L Thompson 
8ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
9ccaff030SJeremy L Thompson // libCEED Operators for MatShell
10ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
11ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator
12ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) {
13ccaff030SJeremy L Thompson   PetscErrorCode ierr;
14ccaff030SJeremy L Thompson   PetscScalar *x, *y;
15d1d35e2fSjeremylt   PetscMemType x_mem_type, y_mem_type;
16ccaff030SJeremy L Thompson 
17ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
18ccaff030SJeremy L Thompson 
19ccaff030SJeremy L Thompson   // Global-to-local
20d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr);
21d1d35e2fSjeremylt   ierr = VecZeroEntries(user->Y_loc); CHKERRQ(ierr);
22ccaff030SJeremy L Thompson 
23ccaff030SJeremy L Thompson   // Setup CEED vectors
24d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
25d1d35e2fSjeremylt                                    &x_mem_type); CHKERRQ(ierr);
26d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr);
27d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
28d1d35e2fSjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
29ccaff030SJeremy L Thompson 
30ccaff030SJeremy L Thompson   // Apply CEED operator
31d1d35e2fSjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed, CEED_REQUEST_IMMEDIATE);
32ccaff030SJeremy L Thompson 
33ccaff030SJeremy L Thompson   // Restore PETSc vectors
34d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
35d1d35e2fSjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
36d1d35e2fSjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
37ccaff030SJeremy L Thompson   CHKERRQ(ierr);
38d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
39ccaff030SJeremy L Thompson 
40ccaff030SJeremy L Thompson   // Local-to-global
41ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
42d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm, user->Y_loc, ADD_VALUES, Y); CHKERRQ(ierr);
43ccaff030SJeremy L Thompson 
44ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
45ccaff030SJeremy L Thompson };
46ccaff030SJeremy L Thompson 
47ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual
48ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
49ccaff030SJeremy L Thompson   PetscErrorCode ierr;
50ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
51ccaff030SJeremy L Thompson 
52ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
53ccaff030SJeremy L Thompson 
54ccaff030SJeremy L Thompson   // Use computed BCs
55d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
56d1d35e2fSjeremylt   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc,
57d1d35e2fSjeremylt                                     user->load_increment, NULL, NULL, NULL);
58ccaff030SJeremy L Thompson   CHKERRQ(ierr);
59ccaff030SJeremy L Thompson 
60ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
61ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
62ccaff030SJeremy L Thompson 
63fe394131Sjeremylt   // Neumann BCs
64d1d35e2fSjeremylt   if (user->neumann_bcs) {
65d1d35e2fSjeremylt     ierr = VecAXPY(Y, -user->load_increment, user->neumann_bcs); CHKERRQ(ierr);
66fe394131Sjeremylt   }
67fe394131Sjeremylt 
68ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
69ccaff030SJeremy L Thompson };
70ccaff030SJeremy L Thompson 
71ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES
72ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
73ccaff030SJeremy L Thompson   PetscErrorCode ierr;
74ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
75ccaff030SJeremy L Thompson 
76ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
77ccaff030SJeremy L Thompson 
7862e9c006SJeremy L Thompson   // Zero boundary values
79d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
80ccaff030SJeremy L Thompson 
81ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
82ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
83ccaff030SJeremy L Thompson 
84ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
85ccaff030SJeremy L Thompson };
86ccaff030SJeremy L Thompson 
87ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian
88ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) {
89ccaff030SJeremy L Thompson   PetscErrorCode ierr;
90ccaff030SJeremy L Thompson   UserMult user;
91ccaff030SJeremy L Thompson 
92ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
93ccaff030SJeremy L Thompson 
94ccaff030SJeremy L Thompson   // Zero boundary values
95ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
96d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
97ccaff030SJeremy L Thompson 
98ccaff030SJeremy L Thompson   // libCEED for local action of Jacobian
99ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
100ccaff030SJeremy L Thompson 
101ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
102ccaff030SJeremy L Thompson };
103ccaff030SJeremy L Thompson 
104ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator
105ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) {
106ccaff030SJeremy L Thompson   PetscErrorCode ierr;
107ccaff030SJeremy L Thompson   UserMultProlongRestr user;
108ccaff030SJeremy L Thompson   PetscScalar *c, *f;
109d1d35e2fSjeremylt   PetscMemType c_mem_type, f_mem_type;
110ccaff030SJeremy L Thompson 
111ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
112ccaff030SJeremy L Thompson 
113ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
114ccaff030SJeremy L Thompson 
115ccaff030SJeremy L Thompson   // Global-to-local
116d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_c); CHKERRQ(ierr);
117d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm_c, X, INSERT_VALUES, user->loc_vec_c);
118ccaff030SJeremy L Thompson   CHKERRQ(ierr);
119d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_f); CHKERRQ(ierr);
120ccaff030SJeremy L Thompson 
121ccaff030SJeremy L Thompson   // Setup CEED vectors
122d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c,
123d1d35e2fSjeremylt                                    &c_mem_type); CHKERRQ(ierr);
124d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->loc_vec_f, &f, &f_mem_type); CHKERRQ(ierr);
125d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER,
126d1d35e2fSjeremylt                      c);
127d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER,
128d1d35e2fSjeremylt                      f);
129ccaff030SJeremy L Thompson 
130ccaff030SJeremy L Thompson   // Apply CEED operator
131d1d35e2fSjeremylt   CeedOperatorApply(user->op_prolong, user->ceed_vec_c, user->ceed_vec_f,
132ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
133ccaff030SJeremy L Thompson 
134ccaff030SJeremy L Thompson   // Restore PETSc vectors
135d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
136d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
137d1d35e2fSjeremylt   ierr = VecRestoreArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c);
138ccaff030SJeremy L Thompson   CHKERRQ(ierr);
139d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->loc_vec_f, &f); CHKERRQ(ierr);
140ccaff030SJeremy L Thompson 
141ccaff030SJeremy L Thompson   // Local-to-global
142ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
143d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm_f, user->loc_vec_f, ADD_VALUES, Y);
144ccaff030SJeremy L Thompson   CHKERRQ(ierr);
145ccaff030SJeremy L Thompson 
146ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
147ccaff030SJeremy L Thompson }
148ccaff030SJeremy L Thompson 
149ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator
150ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
151ccaff030SJeremy L Thompson   PetscErrorCode ierr;
152ccaff030SJeremy L Thompson   UserMultProlongRestr user;
153ccaff030SJeremy L Thompson   PetscScalar *c, *f;
154d1d35e2fSjeremylt   PetscMemType c_mem_type, f_mem_type;
155ccaff030SJeremy L Thompson 
156ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
157ccaff030SJeremy L Thompson 
158ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
159ccaff030SJeremy L Thompson 
160ccaff030SJeremy L Thompson   // Global-to-local
161d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_f); CHKERRQ(ierr);
162d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm_f, X, INSERT_VALUES, user->loc_vec_f);
163ccaff030SJeremy L Thompson   CHKERRQ(ierr);
164d1d35e2fSjeremylt   ierr = VecZeroEntries(user->loc_vec_c); CHKERRQ(ierr);
165ccaff030SJeremy L Thompson 
166ccaff030SJeremy L Thompson   // Setup CEED vectors
167d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f,
168d1d35e2fSjeremylt                                    &f_mem_type); CHKERRQ(ierr);
169d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->loc_vec_c, &c, &c_mem_type); CHKERRQ(ierr);
170d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER,
171d1d35e2fSjeremylt                      f);
172d1d35e2fSjeremylt   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER,
173d1d35e2fSjeremylt                      c);
174ccaff030SJeremy L Thompson 
175ccaff030SJeremy L Thompson   // Apply CEED operator
176d1d35e2fSjeremylt   CeedOperatorApply(user->op_restrict, user->ceed_vec_f, user->ceed_vec_c,
177ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
178ccaff030SJeremy L Thompson 
179ccaff030SJeremy L Thompson   // Restore PETSc vectors
180d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
181d1d35e2fSjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
182d1d35e2fSjeremylt   ierr = VecRestoreArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f);
183ccaff030SJeremy L Thompson   CHKERRQ(ierr);
184d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->loc_vec_c, &c); CHKERRQ(ierr);
185ccaff030SJeremy L Thompson 
186ccaff030SJeremy L Thompson   // Local-to-global
187ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
188d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm_c, user->loc_vec_c, ADD_VALUES, Y);
189ccaff030SJeremy L Thompson   CHKERRQ(ierr);
190ccaff030SJeremy L Thompson 
191ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
192ccaff030SJeremy L Thompson };
1932d93065eSjeremylt 
194ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator
195ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
196ccaff030SJeremy L Thompson   PetscErrorCode ierr;
197ccaff030SJeremy L Thompson   UserMult user;
198ccaff030SJeremy L Thompson 
199ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
200ccaff030SJeremy L Thompson 
201ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
202ccaff030SJeremy L Thompson 
203f7b4142eSJeremy L Thompson   // -- Set physics context
204d1d35e2fSjeremylt   if (user->ctx_phys_smoother)
205d1d35e2fSjeremylt     CeedQFunctionSetContext(user->qf, user->ctx_phys_smoother);
206f7b4142eSJeremy L Thompson 
2072bba3ffaSJeremy L Thompson   // Compute Diagonal via libCEED
2082bba3ffaSJeremy L Thompson   PetscScalar *x;
209d1d35e2fSjeremylt   PetscMemType x_mem_type;
2102bba3ffaSJeremy L Thompson 
2112bba3ffaSJeremy L Thompson   // -- Place PETSc vector in libCEED vector
212d1d35e2fSjeremylt   ierr = VecGetArrayAndMemType(user->X_loc, &x, &x_mem_type); CHKERRQ(ierr);
213d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2142bba3ffaSJeremy L Thompson 
2152bba3ffaSJeremy L Thompson   // -- Compute Diagonal
216d1d35e2fSjeremylt   CeedOperatorLinearAssembleDiagonal(user->op, user->x_ceed,
217ccaff030SJeremy L Thompson                                      CEED_REQUEST_IMMEDIATE);
218ccaff030SJeremy L Thompson 
219f7b4142eSJeremy L Thompson   // -- Reset physics context
220d1d35e2fSjeremylt   if (user->ctx_phys_smoother)
221d1d35e2fSjeremylt     CeedQFunctionSetContext(user->qf, user->ctx_phys);
222f7b4142eSJeremy L Thompson 
223ccaff030SJeremy L Thompson   // -- Local-to-Global
224d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
225d1d35e2fSjeremylt   ierr = VecRestoreArrayAndMemType(user->X_loc, &x); CHKERRQ(ierr);
226ccaff030SJeremy L Thompson   ierr = VecZeroEntries(D); CHKERRQ(ierr);
227d1d35e2fSjeremylt   ierr = DMLocalToGlobal(user->dm, user->X_loc, ADD_VALUES, D); CHKERRQ(ierr);
228ccaff030SJeremy L Thompson 
2292bba3ffaSJeremy L Thompson   // Cleanup
230d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
231ccaff030SJeremy L Thompson 
232ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
233ccaff030SJeremy L Thompson };
2342d93065eSjeremylt 
2352d93065eSjeremylt // This function calculates the strain energy in the final solution
236a3c02c40SJeremy L Thompson PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user,
237d1d35e2fSjeremylt                                    CeedOperator op_energy, Vec X,
238a3c02c40SJeremy L Thompson                                    PetscReal *energy) {
2392d93065eSjeremylt   PetscErrorCode ierr;
2402d93065eSjeremylt   PetscScalar *x;
241d1d35e2fSjeremylt   PetscMemType x_mem_type;
2422d93065eSjeremylt   CeedInt length;
2432d93065eSjeremylt 
2442d93065eSjeremylt   PetscFunctionBeginUser;
2452d93065eSjeremylt 
2462d93065eSjeremylt   // Global-to-local
247d1d35e2fSjeremylt   ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr);
248d1d35e2fSjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr);
249d1d35e2fSjeremylt   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc,
250d1d35e2fSjeremylt                                     user->load_increment, NULL, NULL, NULL);
2515c25879aSJeremy L Thompson   CHKERRQ(ierr);
2522d93065eSjeremylt 
253a3c02c40SJeremy L Thompson   // Setup libCEED input vector
254d1d35e2fSjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
255d1d35e2fSjeremylt                                    &x_mem_type); CHKERRQ(ierr);
256d1d35e2fSjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2572d93065eSjeremylt 
258a3c02c40SJeremy L Thompson   // Setup libCEED output vector
259d1d35e2fSjeremylt   Vec E_loc;
260d1d35e2fSjeremylt   CeedVector e_loc;
261d1d35e2fSjeremylt   ierr = DMCreateLocalVector(dmEnergy, &E_loc); CHKERRQ(ierr);
262d1d35e2fSjeremylt   ierr = VecGetSize(E_loc, &length); CHKERRQ(ierr);
263d1d35e2fSjeremylt   ierr = VecDestroy(&E_loc); CHKERRQ(ierr);
264d1d35e2fSjeremylt   CeedVectorCreate(user->ceed, length, &e_loc);
265a3c02c40SJeremy L Thompson 
2662d93065eSjeremylt   // Apply libCEED operator
267d1d35e2fSjeremylt   CeedOperatorApply(op_energy, user->x_ceed, e_loc, CEED_REQUEST_IMMEDIATE);
2682d93065eSjeremylt 
2692d93065eSjeremylt   // Restore PETSc vector
270d1d35e2fSjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
271d1d35e2fSjeremylt   ierr = VecRestoreArrayRead(user->X_loc, (const PetscScalar **)&x);
2722d93065eSjeremylt   CHKERRQ(ierr);
2732d93065eSjeremylt 
2742d93065eSjeremylt   // Reduce max error
2752d93065eSjeremylt   const CeedScalar *e;
276d1d35e2fSjeremylt   CeedVectorGetArrayRead(e_loc, CEED_MEM_HOST, &e);
2772d93065eSjeremylt   (*energy) = 0;
2782d93065eSjeremylt   for (CeedInt i=0; i<length; i++)
2792d93065eSjeremylt     (*energy) += e[i];
280d1d35e2fSjeremylt   CeedVectorRestoreArrayRead(e_loc, &e);
281d1d35e2fSjeremylt   CeedVectorDestroy(&e_loc);
2822d93065eSjeremylt 
2832d93065eSjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM,
2842d93065eSjeremylt                        user->comm); CHKERRQ(ierr);
2852d93065eSjeremylt 
2862d93065eSjeremylt   PetscFunctionReturn(0);
2872d93065eSjeremylt };
288