xref: /petsc/src/tao/unconstrained/impls/lmvm/tests/ex1.c (revision 5fe01c21d727d5e61b58091b225b37a3f96d9fc7)
1*ec9796c4SHansol Suh const char help[] = "Test TAOLMVM on a least-squares problem";
2*ec9796c4SHansol Suh 
3*ec9796c4SHansol Suh #include <petsctao.h>
4*ec9796c4SHansol Suh #include <petscdevice.h>
5*ec9796c4SHansol Suh 
6*ec9796c4SHansol Suh typedef struct _n_AppCtx {
7*ec9796c4SHansol Suh   Mat A;
8*ec9796c4SHansol Suh   Vec b;
9*ec9796c4SHansol Suh   Vec r;
10*ec9796c4SHansol Suh } AppCtx;
11*ec9796c4SHansol Suh 
LSObjAndGrad(Tao tao,Vec x,PetscReal * obj,Vec g,void * _ctx)12*ec9796c4SHansol Suh static PetscErrorCode LSObjAndGrad(Tao tao, Vec x, PetscReal *obj, Vec g, void *_ctx)
13*ec9796c4SHansol Suh {
14*ec9796c4SHansol Suh   PetscFunctionBegin;
15*ec9796c4SHansol Suh   AppCtx *ctx = (AppCtx *)_ctx;
16*ec9796c4SHansol Suh   PetscCall(VecAXPBY(ctx->r, -1.0, 0.0, ctx->b));
17*ec9796c4SHansol Suh   PetscCall(MatMultAdd(ctx->A, x, ctx->r, ctx->r));
18*ec9796c4SHansol Suh   PetscCall(VecDotRealPart(ctx->r, ctx->r, obj));
19*ec9796c4SHansol Suh   *obj *= 0.5;
20*ec9796c4SHansol Suh   PetscCall(MatMultTranspose(ctx->A, ctx->r, g));
21*ec9796c4SHansol Suh   PetscFunctionReturn(PETSC_SUCCESS);
22*ec9796c4SHansol Suh }
23*ec9796c4SHansol Suh 
main(int argc,char ** argv)24*ec9796c4SHansol Suh int main(int argc, char **argv)
25*ec9796c4SHansol Suh {
26*ec9796c4SHansol Suh   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
27*ec9796c4SHansol Suh   MPI_Comm  comm = PETSC_COMM_WORLD;
28*ec9796c4SHansol Suh   AppCtx    ctx;
29*ec9796c4SHansol Suh   Vec       sol;
30*ec9796c4SHansol Suh   PetscBool flg, cuda = PETSC_FALSE;
31*ec9796c4SHansol Suh 
32*ec9796c4SHansol Suh   PetscInt M = 10;
33*ec9796c4SHansol Suh   PetscInt N = 10;
34*ec9796c4SHansol Suh   PetscOptionsBegin(comm, "", help, "TAO");
35*ec9796c4SHansol Suh   PetscCall(PetscOptionsInt("-m", "data size", NULL, M, &M, NULL));
36*ec9796c4SHansol Suh   PetscCall(PetscOptionsInt("-n", "data size", NULL, N, &N, NULL));
37*ec9796c4SHansol Suh   PetscCall(PetscOptionsGetBool(NULL, NULL, "-cuda", &cuda, &flg));
38*ec9796c4SHansol Suh   PetscOptionsEnd();
39*ec9796c4SHansol Suh 
40*ec9796c4SHansol Suh   if (cuda) {
41*ec9796c4SHansol Suh     VecType vec_type;
42*ec9796c4SHansol Suh     PetscCall(VecCreateSeqCUDA(comm, N, &ctx.b));
43*ec9796c4SHansol Suh     PetscCall(VecGetType(ctx.b, &vec_type));
44*ec9796c4SHansol Suh     PetscCall(MatCreateDenseFromVecType(comm, vec_type, M, N, PETSC_DECIDE, PETSC_DECIDE, -1, NULL, &ctx.A));
45*ec9796c4SHansol Suh     PetscCall(MatCreateVecs(ctx.A, &sol, NULL));
46*ec9796c4SHansol Suh   } else {
47*ec9796c4SHansol Suh     PetscCall(MatCreateDense(comm, PETSC_DECIDE, PETSC_DECIDE, M, N, NULL, &ctx.A));
48*ec9796c4SHansol Suh     PetscCall(MatCreateVecs(ctx.A, &sol, &ctx.b));
49*ec9796c4SHansol Suh   }
50*ec9796c4SHansol Suh   PetscCall(VecDuplicate(ctx.b, &ctx.r));
51*ec9796c4SHansol Suh   PetscCall(VecZeroEntries(sol));
52*ec9796c4SHansol Suh 
53*ec9796c4SHansol Suh   PetscRandom rand;
54*ec9796c4SHansol Suh   PetscCall(PetscRandomCreate(comm, &rand));
55*ec9796c4SHansol Suh   PetscCall(PetscRandomSetFromOptions(rand));
56*ec9796c4SHansol Suh   PetscCall(MatSetRandom(ctx.A, rand));
57*ec9796c4SHansol Suh   PetscCall(VecSetRandom(ctx.b, rand));
58*ec9796c4SHansol Suh   PetscCall(PetscRandomDestroy(&rand));
59*ec9796c4SHansol Suh 
60*ec9796c4SHansol Suh   Tao tao;
61*ec9796c4SHansol Suh   PetscCall(TaoCreate(comm, &tao));
62*ec9796c4SHansol Suh   PetscCall(TaoSetSolution(tao, sol));
63*ec9796c4SHansol Suh   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, LSObjAndGrad, &ctx));
64*ec9796c4SHansol Suh   PetscCall(TaoSetType(tao, TAOLMVM));
65*ec9796c4SHansol Suh   PetscCall(TaoSetFromOptions(tao));
66*ec9796c4SHansol Suh   PetscCall(TaoSolve(tao));
67*ec9796c4SHansol Suh   PetscCall(TaoDestroy(&tao));
68*ec9796c4SHansol Suh 
69*ec9796c4SHansol Suh   PetscCall(VecDestroy(&ctx.r));
70*ec9796c4SHansol Suh   PetscCall(VecDestroy(&sol));
71*ec9796c4SHansol Suh   PetscCall(VecDestroy(&ctx.b));
72*ec9796c4SHansol Suh   PetscCall(MatDestroy(&ctx.A));
73*ec9796c4SHansol Suh 
74*ec9796c4SHansol Suh   PetscCall(PetscFinalize());
75*ec9796c4SHansol Suh   return 0;
76*ec9796c4SHansol Suh }
77*ec9796c4SHansol Suh 
78*ec9796c4SHansol Suh /*TEST
79*ec9796c4SHansol Suh 
80*ec9796c4SHansol Suh   build:
81*ec9796c4SHansol Suh     requires: !complex !__float128 !single !defined(PETSC_USE_64BIT_INDICES)
82*ec9796c4SHansol Suh 
83*ec9796c4SHansol Suh   test:
84*ec9796c4SHansol Suh     suffix: 0
85*ec9796c4SHansol Suh     args: -tao_monitor -tao_ls_gtol 1.e-6 -tao_view -tao_lmvm_mat_lmvm_hist_size 20 -tao_ls_type more-thuente -tao_lmvm_mat_lmvm_scale_type none -tao_lmvm_mat_type lmvmbfgs
86*ec9796c4SHansol Suh 
87*ec9796c4SHansol Suh   test:
88*ec9796c4SHansol Suh     suffix: 1
89*ec9796c4SHansol Suh     args: -tao_monitor -tao_ls_gtol 1.e-6 -tao_view -tao_lmvm_mat_lmvm_hist_size 20 -tao_ls_type more-thuente -tao_lmvm_mat_lmvm_scale_type none -tao_lmvm_mat_type lmvmdbfgs
90*ec9796c4SHansol Suh 
91*ec9796c4SHansol Suh   test:
92*ec9796c4SHansol Suh     suffix: 2
93*ec9796c4SHansol Suh     args: -tao_monitor -tao_ls_gtol 1.e-6 -tao_view -tao_lmvm_mat_lmvm_hist_size 20 -tao_ls_type more-thuente -tao_lmvm_mat_type lmvmdbfgs -tao_lmvm_mat_lmvm_scale_type none -tao_lmvm_mat_lbfgs_type {{inplace reorder}}
94*ec9796c4SHansol Suh 
95*ec9796c4SHansol Suh TEST*/
96