xref: /libCEED/examples/petsc/src/matops.c (revision f1c40530bc6f02bc6ca847af718b8b41cc3fd3ca)
1e83e87a5Sjeremylt #include "../include/matops.h"
2e83e87a5Sjeremylt #include "../include/petscutils.h"
3e83e87a5Sjeremylt 
4e83e87a5Sjeremylt // -----------------------------------------------------------------------------
5e83e87a5Sjeremylt // This function returns the computed diagonal of the operator
6e83e87a5Sjeremylt // -----------------------------------------------------------------------------
7e83e87a5Sjeremylt PetscErrorCode MatGetDiag(Mat A, Vec D) {
8e83e87a5Sjeremylt   PetscErrorCode ierr;
9e83e87a5Sjeremylt   UserO user;
10e83e87a5Sjeremylt 
11e83e87a5Sjeremylt   PetscFunctionBeginUser;
12e83e87a5Sjeremylt 
13e83e87a5Sjeremylt   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
14e83e87a5Sjeremylt 
15e83e87a5Sjeremylt   // Compute Diagonal via libCEED
16*f1c40530SJeremy L Thompson   PetscScalar *y;
179b072555Sjeremylt   PetscMemType mem_type;
18e83e87a5Sjeremylt 
19e83e87a5Sjeremylt   // -- Place PETSc vector in libCEED vector
20*f1c40530SJeremy L Thompson   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &mem_type); CHKERRQ(ierr);
21*f1c40530SJeremy L Thompson   CeedVectorSetArray(user->y_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, y);
22e83e87a5Sjeremylt 
23e83e87a5Sjeremylt   // -- Compute Diagonal
24*f1c40530SJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op, user->y_ceed,
25e83e87a5Sjeremylt                                      CEED_REQUEST_IMMEDIATE);
26e83e87a5Sjeremylt 
27e83e87a5Sjeremylt   // -- Local-to-Global
28*f1c40530SJeremy L Thompson   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(mem_type), NULL);
29*f1c40530SJeremy L Thompson   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
30e83e87a5Sjeremylt   ierr = VecZeroEntries(D); CHKERRQ(ierr);
31*f1c40530SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, user->Y_loc, ADD_VALUES, D); CHKERRQ(ierr);
32e83e87a5Sjeremylt 
33e83e87a5Sjeremylt   PetscFunctionReturn(0);
34e83e87a5Sjeremylt };
35e83e87a5Sjeremylt 
36e83e87a5Sjeremylt // -----------------------------------------------------------------------------
37e83e87a5Sjeremylt // This function uses libCEED to compute the action of the Laplacian with
38e83e87a5Sjeremylt // Dirichlet boundary conditions
39e83e87a5Sjeremylt // -----------------------------------------------------------------------------
40e83e87a5Sjeremylt PetscErrorCode ApplyLocal_Ceed(Vec X, Vec Y, UserO user) {
41e83e87a5Sjeremylt   PetscErrorCode ierr;
42e83e87a5Sjeremylt   PetscScalar *x, *y;
439b072555Sjeremylt   PetscMemType x_mem_type, y_mem_type;
44e83e87a5Sjeremylt 
45e83e87a5Sjeremylt   PetscFunctionBeginUser;
46e83e87a5Sjeremylt 
47e83e87a5Sjeremylt   // Global-to-local
489b072555Sjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr);
49e83e87a5Sjeremylt 
50e83e87a5Sjeremylt   // Setup libCEED vectors
519b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
529b072555Sjeremylt                                    &x_mem_type); CHKERRQ(ierr);
539b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr);
549b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
559b072555Sjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
56e83e87a5Sjeremylt 
57e83e87a5Sjeremylt   // Apply libCEED operator
589b072555Sjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed, CEED_REQUEST_IMMEDIATE);
59e83e87a5Sjeremylt 
60e83e87a5Sjeremylt   // Restore PETSc vectors
619b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
629b072555Sjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
639b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
64e83e87a5Sjeremylt   CHKERRQ(ierr);
659b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
66e83e87a5Sjeremylt 
67e83e87a5Sjeremylt   // Local-to-global
68e83e87a5Sjeremylt   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
699b072555Sjeremylt   ierr = DMLocalToGlobal(user->dm, user->Y_loc, ADD_VALUES, Y); CHKERRQ(ierr);
70e83e87a5Sjeremylt 
71e83e87a5Sjeremylt   PetscFunctionReturn(0);
72e83e87a5Sjeremylt };
73e83e87a5Sjeremylt 
74e83e87a5Sjeremylt // -----------------------------------------------------------------------------
75e83e87a5Sjeremylt // This function wraps the libCEED operator for a MatShell
76e83e87a5Sjeremylt // -----------------------------------------------------------------------------
77e83e87a5Sjeremylt PetscErrorCode MatMult_Ceed(Mat A, Vec X, Vec Y) {
78e83e87a5Sjeremylt   PetscErrorCode ierr;
79e83e87a5Sjeremylt   UserO user;
80e83e87a5Sjeremylt 
81e83e87a5Sjeremylt   PetscFunctionBeginUser;
82e83e87a5Sjeremylt 
83e83e87a5Sjeremylt   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
84e83e87a5Sjeremylt 
85e83e87a5Sjeremylt   // libCEED for local action of residual evaluator
86e83e87a5Sjeremylt   ierr = ApplyLocal_Ceed(X, Y, user); CHKERRQ(ierr);
87e83e87a5Sjeremylt 
88e83e87a5Sjeremylt   PetscFunctionReturn(0);
89e83e87a5Sjeremylt };
90e83e87a5Sjeremylt 
91e83e87a5Sjeremylt // -----------------------------------------------------------------------------
92e83e87a5Sjeremylt // This function wraps the libCEED operator for a SNES residual evaluation
93e83e87a5Sjeremylt // -----------------------------------------------------------------------------
94e83e87a5Sjeremylt PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
95e83e87a5Sjeremylt   PetscErrorCode ierr;
96e83e87a5Sjeremylt   UserO user = (UserO)ctx;
97e83e87a5Sjeremylt 
98e83e87a5Sjeremylt   PetscFunctionBeginUser;
99e83e87a5Sjeremylt 
100e83e87a5Sjeremylt   // libCEED for local action of residual evaluator
101e83e87a5Sjeremylt   ierr = ApplyLocal_Ceed(X, Y, user); CHKERRQ(ierr);
102e83e87a5Sjeremylt 
103e83e87a5Sjeremylt   PetscFunctionReturn(0);
104e83e87a5Sjeremylt };
105e83e87a5Sjeremylt 
106e83e87a5Sjeremylt // -----------------------------------------------------------------------------
107e83e87a5Sjeremylt // This function uses libCEED to compute the action of the prolongation operator
108e83e87a5Sjeremylt // -----------------------------------------------------------------------------
109e83e87a5Sjeremylt PetscErrorCode MatMult_Prolong(Mat A, Vec X, Vec Y) {
110e83e87a5Sjeremylt   PetscErrorCode ierr;
111e83e87a5Sjeremylt   UserProlongRestr user;
112e83e87a5Sjeremylt   PetscScalar *c, *f;
1139b072555Sjeremylt   PetscMemType c_mem_type, f_mem_type;
114e83e87a5Sjeremylt 
115e83e87a5Sjeremylt   PetscFunctionBeginUser;
116e83e87a5Sjeremylt 
117e83e87a5Sjeremylt   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
118e83e87a5Sjeremylt 
119e83e87a5Sjeremylt   // Global-to-local
1209b072555Sjeremylt   ierr = VecZeroEntries(user->loc_vec_c); CHKERRQ(ierr);
1219b072555Sjeremylt   ierr = DMGlobalToLocal(user->dmc, X, INSERT_VALUES, user->loc_vec_c);
122e83e87a5Sjeremylt   CHKERRQ(ierr);
123e83e87a5Sjeremylt 
124e83e87a5Sjeremylt   // Setup libCEED vectors
1259b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c,
1269b072555Sjeremylt                                    &c_mem_type); CHKERRQ(ierr);
1279b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->loc_vec_f, &f, &f_mem_type); CHKERRQ(ierr);
1289b072555Sjeremylt   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER,
1299b072555Sjeremylt                      c);
1309b072555Sjeremylt   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER,
1319b072555Sjeremylt                      f);
132e83e87a5Sjeremylt 
133e83e87a5Sjeremylt   // Apply libCEED operator
1349b072555Sjeremylt   CeedOperatorApply(user->op_prolong, user->ceed_vec_c, user->ceed_vec_f,
135e83e87a5Sjeremylt                     CEED_REQUEST_IMMEDIATE);
136e83e87a5Sjeremylt 
137e83e87a5Sjeremylt   // Restore PETSc vectors
1389b072555Sjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
1399b072555Sjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
1409b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c);
141e83e87a5Sjeremylt   CHKERRQ(ierr);
1429b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(user->loc_vec_f, &f); CHKERRQ(ierr);
143e83e87a5Sjeremylt 
144e83e87a5Sjeremylt   // Multiplicity
1459b072555Sjeremylt   ierr = VecPointwiseMult(user->loc_vec_f, user->loc_vec_f, user->mult_vec);
146e83e87a5Sjeremylt 
147e83e87a5Sjeremylt   // Local-to-global
148e83e87a5Sjeremylt   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
1499b072555Sjeremylt   ierr = DMLocalToGlobal(user->dmf, user->loc_vec_f, ADD_VALUES, Y);
150e83e87a5Sjeremylt   CHKERRQ(ierr);
151e83e87a5Sjeremylt 
152e83e87a5Sjeremylt   PetscFunctionReturn(0);
153e83e87a5Sjeremylt };
154e83e87a5Sjeremylt 
155e83e87a5Sjeremylt // -----------------------------------------------------------------------------
156e83e87a5Sjeremylt // This function uses libCEED to compute the action of the restriction operator
157e83e87a5Sjeremylt // -----------------------------------------------------------------------------
158e83e87a5Sjeremylt PetscErrorCode MatMult_Restrict(Mat A, Vec X, Vec Y) {
159e83e87a5Sjeremylt   PetscErrorCode ierr;
160e83e87a5Sjeremylt   UserProlongRestr user;
161e83e87a5Sjeremylt   PetscScalar *c, *f;
1629b072555Sjeremylt   PetscMemType c_mem_type, f_mem_type;
163e83e87a5Sjeremylt 
164e83e87a5Sjeremylt   PetscFunctionBeginUser;
165e83e87a5Sjeremylt 
166e83e87a5Sjeremylt   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
167e83e87a5Sjeremylt 
168e83e87a5Sjeremylt   // Global-to-local
1699b072555Sjeremylt   ierr = VecZeroEntries(user->loc_vec_f); CHKERRQ(ierr);
1709b072555Sjeremylt   ierr = DMGlobalToLocal(user->dmf, X, INSERT_VALUES, user->loc_vec_f);
171e83e87a5Sjeremylt   CHKERRQ(ierr);
172e83e87a5Sjeremylt 
173e83e87a5Sjeremylt   // Multiplicity
1749b072555Sjeremylt   ierr = VecPointwiseMult(user->loc_vec_f, user->loc_vec_f, user->mult_vec);
175e83e87a5Sjeremylt   CHKERRQ(ierr);
176e83e87a5Sjeremylt 
177e83e87a5Sjeremylt   // Setup libCEED vectors
1789b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f,
1799b072555Sjeremylt                                    &f_mem_type); CHKERRQ(ierr);
1809b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->loc_vec_c, &c, &c_mem_type); CHKERRQ(ierr);
1819b072555Sjeremylt   CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER,
1829b072555Sjeremylt                      f);
1839b072555Sjeremylt   CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER,
1849b072555Sjeremylt                      c);
185e83e87a5Sjeremylt 
186e83e87a5Sjeremylt   // Apply CEED operator
1879b072555Sjeremylt   CeedOperatorApply(user->op_restrict, user->ceed_vec_f, user->ceed_vec_c,
188e83e87a5Sjeremylt                     CEED_REQUEST_IMMEDIATE);
189e83e87a5Sjeremylt 
190e83e87a5Sjeremylt   // Restore PETSc vectors
1919b072555Sjeremylt   CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL);
1929b072555Sjeremylt   CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL);
1939b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f);
194e83e87a5Sjeremylt   CHKERRQ(ierr);
1959b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(user->loc_vec_c, &c); CHKERRQ(ierr);
196e83e87a5Sjeremylt 
197e83e87a5Sjeremylt   // Local-to-global
198e83e87a5Sjeremylt   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
1999b072555Sjeremylt   ierr = DMLocalToGlobal(user->dmc, user->loc_vec_c, ADD_VALUES, Y);
200e83e87a5Sjeremylt   CHKERRQ(ierr);
201e83e87a5Sjeremylt 
202e83e87a5Sjeremylt   PetscFunctionReturn(0);
203e83e87a5Sjeremylt };
204e83e87a5Sjeremylt 
205e83e87a5Sjeremylt // -----------------------------------------------------------------------------
206e83e87a5Sjeremylt // This function calculates the error in the final solution
207e83e87a5Sjeremylt // -----------------------------------------------------------------------------
2089b072555Sjeremylt PetscErrorCode ComputeErrorMax(UserO user, CeedOperator op_error,
209e83e87a5Sjeremylt                                Vec X, CeedVector target,
2109b072555Sjeremylt                                PetscScalar *max_error) {
211e83e87a5Sjeremylt   PetscErrorCode ierr;
212e83e87a5Sjeremylt   PetscScalar *x;
2139b072555Sjeremylt   PetscMemType mem_type;
214e83e87a5Sjeremylt   CeedVector collocated_error;
2151f9221feSJeremy L Thompson   CeedSize length;
216e83e87a5Sjeremylt 
217e83e87a5Sjeremylt   PetscFunctionBeginUser;
218e83e87a5Sjeremylt   CeedVectorGetLength(target, &length);
219e83e87a5Sjeremylt   CeedVectorCreate(user->ceed, length, &collocated_error);
220e83e87a5Sjeremylt 
221e83e87a5Sjeremylt   // Global-to-local
2229b072555Sjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr);
223e83e87a5Sjeremylt 
224e83e87a5Sjeremylt   // Setup CEED vector
2259b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->X_loc, &x, &mem_type); CHKERRQ(ierr);
2269b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x);
227e83e87a5Sjeremylt 
228e83e87a5Sjeremylt   // Apply CEED operator
2299b072555Sjeremylt   CeedOperatorApply(op_error, user->x_ceed, collocated_error,
230e83e87a5Sjeremylt                     CEED_REQUEST_IMMEDIATE);
231e83e87a5Sjeremylt 
232e83e87a5Sjeremylt   // Restore PETSc vector
2339b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(mem_type), NULL);
2349b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
235e83e87a5Sjeremylt   CHKERRQ(ierr);
236e83e87a5Sjeremylt 
237e83e87a5Sjeremylt   // Reduce max error
2389b072555Sjeremylt   *max_error = 0;
239e83e87a5Sjeremylt   const CeedScalar *e;
240e83e87a5Sjeremylt   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
241e83e87a5Sjeremylt   for (CeedInt i=0; i<length; i++) {
2429b072555Sjeremylt     *max_error = PetscMax(*max_error, PetscAbsScalar(e[i]));
243e83e87a5Sjeremylt   }
244e83e87a5Sjeremylt   CeedVectorRestoreArrayRead(collocated_error, &e);
2459b072555Sjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX,
246e83e87a5Sjeremylt                        user->comm);
247e83e87a5Sjeremylt   CHKERRQ(ierr);
248e83e87a5Sjeremylt 
249e83e87a5Sjeremylt   // Cleanup
250e83e87a5Sjeremylt   CeedVectorDestroy(&collocated_error);
251e83e87a5Sjeremylt 
252e83e87a5Sjeremylt   PetscFunctionReturn(0);
253e83e87a5Sjeremylt };
254e83e87a5Sjeremylt 
255e83e87a5Sjeremylt // -----------------------------------------------------------------------------
256