1! Program usage: mpiexec -n 1 chwirut1f [-help] [all TAO options] 2! 3! Description: This example demonstrates use of the TAO package to solve a 4! nonlinear least-squares problem on a single processor. We minimize the 5! Chwirut function: 6! sum_{i=0}^{n/2-1} ( alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) 7! 8! The C version of this code is chwirut1.c 9! 10#include <petsc/finclude/petsctao.h> 11module chwirut2fmodule 12 use petscmpi ! or mpi or mpi_f08 13 use petsctao 14 implicit none 15 PetscReal t(0:213) 16 PetscReal y(0:213) 17 PetscInt, parameter :: m = 214, n = 3 18 PetscMPIInt, parameter :: nn = n 19 PetscMPIInt rank 20 PetscMPIInt size 21 PetscMPIInt, parameter :: idle_tag = 2000, die_tag = 3000 22 PetscMPIInt, parameter :: zero = 0, one = 1 23 24contains 25 subroutine InitializeData() 26 27 PetscInt i 28 i = 0 29 y(i) = 92.9000; t(i) = 0.5000; i = i + 1 30 y(i) = 78.7000; t(i) = 0.6250; i = i + 1 31 y(i) = 64.2000; t(i) = 0.7500; i = i + 1 32 y(i) = 64.9000; t(i) = 0.8750; i = i + 1 33 y(i) = 57.1000; t(i) = 1.0000; i = i + 1 34 y(i) = 43.3000; t(i) = 1.2500; i = i + 1 35 y(i) = 31.1000; t(i) = 1.7500; i = i + 1 36 y(i) = 23.6000; t(i) = 2.2500; i = i + 1 37 y(i) = 31.0500; t(i) = 1.7500; i = i + 1 38 y(i) = 23.7750; t(i) = 2.2500; i = i + 1 39 y(i) = 17.7375; t(i) = 2.7500; i = i + 1 40 y(i) = 13.8000; t(i) = 3.2500; i = i + 1 41 y(i) = 11.5875; t(i) = 3.7500; i = i + 1 42 y(i) = 9.4125; t(i) = 4.2500; i = i + 1 43 y(i) = 7.7250; t(i) = 4.7500; i = i + 1 44 y(i) = 7.3500; t(i) = 5.2500; i = i + 1 45 y(i) = 8.0250; t(i) = 5.7500; i = i + 1 46 y(i) = 90.6000; t(i) = 0.5000; i = i + 1 47 y(i) = 76.9000; t(i) = 0.6250; i = i + 1 48 y(i) = 71.6000; t(i) = 0.7500; i = i + 1 49 y(i) = 63.6000; t(i) = 0.8750; i = i + 1 50 y(i) = 54.0000; t(i) = 1.0000; i = i + 1 51 y(i) = 39.2000; t(i) = 1.2500; i = i + 1 52 y(i) = 29.3000; t(i) = 1.7500; i = i + 1 53 y(i) = 21.4000; t(i) = 2.2500; i = i + 1 54 y(i) = 29.1750; t(i) = 1.7500; i = i + 1 55 y(i) = 22.1250; t(i) = 2.2500; i = i + 1 56 y(i) = 17.5125; t(i) = 2.7500; i = i + 1 57 y(i) = 14.2500; t(i) = 3.2500; i = i + 1 58 y(i) = 9.4500; t(i) = 3.7500; i = i + 1 59 y(i) = 9.1500; t(i) = 4.2500; i = i + 1 60 y(i) = 7.9125; t(i) = 4.7500; i = i + 1 61 y(i) = 8.4750; t(i) = 5.2500; i = i + 1 62 y(i) = 6.1125; t(i) = 5.7500; i = i + 1 63 y(i) = 80.0000; t(i) = 0.5000; i = i + 1 64 y(i) = 79.0000; t(i) = 0.6250; i = i + 1 65 y(i) = 63.8000; t(i) = 0.7500; i = i + 1 66 y(i) = 57.2000; t(i) = 0.8750; i = i + 1 67 y(i) = 53.2000; t(i) = 1.0000; i = i + 1 68 y(i) = 42.5000; t(i) = 1.2500; i = i + 1 69 y(i) = 26.8000; t(i) = 1.7500; i = i + 1 70 y(i) = 20.4000; t(i) = 2.2500; i = i + 1 71 y(i) = 26.8500; t(i) = 1.7500; i = i + 1 72 y(i) = 21.0000; t(i) = 2.2500; i = i + 1 73 y(i) = 16.4625; t(i) = 2.7500; i = i + 1 74 y(i) = 12.5250; t(i) = 3.2500; i = i + 1 75 y(i) = 10.5375; t(i) = 3.7500; i = i + 1 76 y(i) = 8.5875; t(i) = 4.2500; i = i + 1 77 y(i) = 7.1250; t(i) = 4.7500; i = i + 1 78 y(i) = 6.1125; t(i) = 5.2500; i = i + 1 79 y(i) = 5.9625; t(i) = 5.7500; i = i + 1 80 y(i) = 74.1000; t(i) = 0.5000; i = i + 1 81 y(i) = 67.3000; t(i) = 0.6250; i = i + 1 82 y(i) = 60.8000; t(i) = 0.7500; i = i + 1 83 y(i) = 55.5000; t(i) = 0.8750; i = i + 1 84 y(i) = 50.3000; t(i) = 1.0000; i = i + 1 85 y(i) = 41.0000; t(i) = 1.2500; i = i + 1 86 y(i) = 29.4000; t(i) = 1.7500; i = i + 1 87 y(i) = 20.4000; t(i) = 2.2500; i = i + 1 88 y(i) = 29.3625; t(i) = 1.7500; i = i + 1 89 y(i) = 21.1500; t(i) = 2.2500; i = i + 1 90 y(i) = 16.7625; t(i) = 2.7500; i = i + 1 91 y(i) = 13.2000; t(i) = 3.2500; i = i + 1 92 y(i) = 10.8750; t(i) = 3.7500; i = i + 1 93 y(i) = 8.1750; t(i) = 4.2500; i = i + 1 94 y(i) = 7.3500; t(i) = 4.7500; i = i + 1 95 y(i) = 5.9625; t(i) = 5.2500; i = i + 1 96 y(i) = 5.6250; t(i) = 5.7500; i = i + 1 97 y(i) = 81.5000; t(i) = .5000; i = i + 1 98 y(i) = 62.4000; t(i) = .7500; i = i + 1 99 y(i) = 32.5000; t(i) = 1.5000; i = i + 1 100 y(i) = 12.4100; t(i) = 3.0000; i = i + 1 101 y(i) = 13.1200; t(i) = 3.0000; i = i + 1 102 y(i) = 15.5600; t(i) = 3.0000; i = i + 1 103 y(i) = 5.6300; t(i) = 6.0000; i = i + 1 104 y(i) = 78.0000; t(i) = .5000; i = i + 1 105 y(i) = 59.9000; t(i) = .7500; i = i + 1 106 y(i) = 33.2000; t(i) = 1.5000; i = i + 1 107 y(i) = 13.8400; t(i) = 3.0000; i = i + 1 108 y(i) = 12.7500; t(i) = 3.0000; i = i + 1 109 y(i) = 14.6200; t(i) = 3.0000; i = i + 1 110 y(i) = 3.9400; t(i) = 6.0000; i = i + 1 111 y(i) = 76.8000; t(i) = .5000; i = i + 1 112 y(i) = 61.0000; t(i) = .7500; i = i + 1 113 y(i) = 32.9000; t(i) = 1.5000; i = i + 1 114 y(i) = 13.8700; t(i) = 3.0000; i = i + 1 115 y(i) = 11.8100; t(i) = 3.0000; i = i + 1 116 y(i) = 13.3100; t(i) = 3.0000; i = i + 1 117 y(i) = 5.4400; t(i) = 6.0000; i = i + 1 118 y(i) = 78.0000; t(i) = .5000; i = i + 1 119 y(i) = 63.5000; t(i) = .7500; i = i + 1 120 y(i) = 33.8000; t(i) = 1.5000; i = i + 1 121 y(i) = 12.5600; t(i) = 3.0000; i = i + 1 122 y(i) = 5.6300; t(i) = 6.0000; i = i + 1 123 y(i) = 12.7500; t(i) = 3.0000; i = i + 1 124 y(i) = 13.1200; t(i) = 3.0000; i = i + 1 125 y(i) = 5.4400; t(i) = 6.0000; i = i + 1 126 y(i) = 76.8000; t(i) = .5000; i = i + 1 127 y(i) = 60.0000; t(i) = .7500; i = i + 1 128 y(i) = 47.8000; t(i) = 1.0000; i = i + 1 129 y(i) = 32.0000; t(i) = 1.5000; i = i + 1 130 y(i) = 22.2000; t(i) = 2.0000; i = i + 1 131 y(i) = 22.5700; t(i) = 2.0000; i = i + 1 132 y(i) = 18.8200; t(i) = 2.5000; i = i + 1 133 y(i) = 13.9500; t(i) = 3.0000; i = i + 1 134 y(i) = 11.2500; t(i) = 4.0000; i = i + 1 135 y(i) = 9.0000; t(i) = 5.0000; i = i + 1 136 y(i) = 6.6700; t(i) = 6.0000; i = i + 1 137 y(i) = 75.8000; t(i) = .5000; i = i + 1 138 y(i) = 62.0000; t(i) = .7500; i = i + 1 139 y(i) = 48.8000; t(i) = 1.0000; i = i + 1 140 y(i) = 35.2000; t(i) = 1.5000; i = i + 1 141 y(i) = 20.0000; t(i) = 2.0000; i = i + 1 142 y(i) = 20.3200; t(i) = 2.0000; i = i + 1 143 y(i) = 19.3100; t(i) = 2.5000; i = i + 1 144 y(i) = 12.7500; t(i) = 3.0000; i = i + 1 145 y(i) = 10.4200; t(i) = 4.0000; i = i + 1 146 y(i) = 7.3100; t(i) = 5.0000; i = i + 1 147 y(i) = 7.4200; t(i) = 6.0000; i = i + 1 148 y(i) = 70.5000; t(i) = .5000; i = i + 1 149 y(i) = 59.5000; t(i) = .7500; i = i + 1 150 y(i) = 48.5000; t(i) = 1.0000; i = i + 1 151 y(i) = 35.8000; t(i) = 1.5000; i = i + 1 152 y(i) = 21.0000; t(i) = 2.0000; i = i + 1 153 y(i) = 21.6700; t(i) = 2.0000; i = i + 1 154 y(i) = 21.0000; t(i) = 2.5000; i = i + 1 155 y(i) = 15.6400; t(i) = 3.0000; i = i + 1 156 y(i) = 8.1700; t(i) = 4.0000; i = i + 1 157 y(i) = 8.5500; t(i) = 5.0000; i = i + 1 158 y(i) = 10.1200; t(i) = 6.0000; i = i + 1 159 y(i) = 78.0000; t(i) = .5000; i = i + 1 160 y(i) = 66.0000; t(i) = .6250; i = i + 1 161 y(i) = 62.0000; t(i) = .7500; i = i + 1 162 y(i) = 58.0000; t(i) = .8750; i = i + 1 163 y(i) = 47.7000; t(i) = 1.0000; i = i + 1 164 y(i) = 37.8000; t(i) = 1.2500; i = i + 1 165 y(i) = 20.2000; t(i) = 2.2500; i = i + 1 166 y(i) = 21.0700; t(i) = 2.2500; i = i + 1 167 y(i) = 13.8700; t(i) = 2.7500; i = i + 1 168 y(i) = 9.6700; t(i) = 3.2500; i = i + 1 169 y(i) = 7.7600; t(i) = 3.7500; i = i + 1 170 y(i) = 5.4400; t(i) = 4.2500; i = i + 1 171 y(i) = 4.8700; t(i) = 4.7500; i = i + 1 172 y(i) = 4.0100; t(i) = 5.2500; i = i + 1 173 y(i) = 3.7500; t(i) = 5.7500; i = i + 1 174 y(i) = 24.1900; t(i) = 3.0000; i = i + 1 175 y(i) = 25.7600; t(i) = 3.0000; i = i + 1 176 y(i) = 18.0700; t(i) = 3.0000; i = i + 1 177 y(i) = 11.8100; t(i) = 3.0000; i = i + 1 178 y(i) = 12.0700; t(i) = 3.0000; i = i + 1 179 y(i) = 16.1200; t(i) = 3.0000; i = i + 1 180 y(i) = 70.8000; t(i) = .5000; i = i + 1 181 y(i) = 54.7000; t(i) = .7500; i = i + 1 182 y(i) = 48.0000; t(i) = 1.0000; i = i + 1 183 y(i) = 39.8000; t(i) = 1.5000; i = i + 1 184 y(i) = 29.8000; t(i) = 2.0000; i = i + 1 185 y(i) = 23.7000; t(i) = 2.5000; i = i + 1 186 y(i) = 29.6200; t(i) = 2.0000; i = i + 1 187 y(i) = 23.8100; t(i) = 2.5000; i = i + 1 188 y(i) = 17.7000; t(i) = 3.0000; i = i + 1 189 y(i) = 11.5500; t(i) = 4.0000; i = i + 1 190 y(i) = 12.0700; t(i) = 5.0000; i = i + 1 191 y(i) = 8.7400; t(i) = 6.0000; i = i + 1 192 y(i) = 80.7000; t(i) = .5000; i = i + 1 193 y(i) = 61.3000; t(i) = .7500; i = i + 1 194 y(i) = 47.5000; t(i) = 1.0000; i = i + 1 195 y(i) = 29.0000; t(i) = 1.5000; i = i + 1 196 y(i) = 24.0000; t(i) = 2.0000; i = i + 1 197 y(i) = 17.7000; t(i) = 2.5000; i = i + 1 198 y(i) = 24.5600; t(i) = 2.0000; i = i + 1 199 y(i) = 18.6700; t(i) = 2.5000; i = i + 1 200 y(i) = 16.2400; t(i) = 3.0000; i = i + 1 201 y(i) = 8.7400; t(i) = 4.0000; i = i + 1 202 y(i) = 7.8700; t(i) = 5.0000; i = i + 1 203 y(i) = 8.5100; t(i) = 6.0000; i = i + 1 204 y(i) = 66.7000; t(i) = .5000; i = i + 1 205 y(i) = 59.2000; t(i) = .7500; i = i + 1 206 y(i) = 40.8000; t(i) = 1.0000; i = i + 1 207 y(i) = 30.7000; t(i) = 1.5000; i = i + 1 208 y(i) = 25.7000; t(i) = 2.0000; i = i + 1 209 y(i) = 16.3000; t(i) = 2.5000; i = i + 1 210 y(i) = 25.9900; t(i) = 2.0000; i = i + 1 211 y(i) = 16.9500; t(i) = 2.5000; i = i + 1 212 y(i) = 13.3500; t(i) = 3.0000; i = i + 1 213 y(i) = 8.6200; t(i) = 4.0000; i = i + 1 214 y(i) = 7.2000; t(i) = 5.0000; i = i + 1 215 y(i) = 6.6400; t(i) = 6.0000; i = i + 1 216 y(i) = 13.6900; t(i) = 3.0000; i = i + 1 217 y(i) = 81.0000; t(i) = .5000; i = i + 1 218 y(i) = 64.5000; t(i) = .7500; i = i + 1 219 y(i) = 35.5000; t(i) = 1.5000; i = i + 1 220 y(i) = 13.3100; t(i) = 3.0000; i = i + 1 221 y(i) = 4.8700; t(i) = 6.0000; i = i + 1 222 y(i) = 12.9400; t(i) = 3.0000; i = i + 1 223 y(i) = 5.0600; t(i) = 6.0000; i = i + 1 224 y(i) = 15.1900; t(i) = 3.0000; i = i + 1 225 y(i) = 14.6200; t(i) = 3.0000; i = i + 1 226 y(i) = 15.6400; t(i) = 3.0000; i = i + 1 227 y(i) = 25.5000; t(i) = 1.7500; i = i + 1 228 y(i) = 25.9500; t(i) = 1.7500; i = i + 1 229 y(i) = 81.7000; t(i) = .5000; i = i + 1 230 y(i) = 61.6000; t(i) = .7500; i = i + 1 231 y(i) = 29.8000; t(i) = 1.7500; i = i + 1 232 y(i) = 29.8100; t(i) = 1.7500; i = i + 1 233 y(i) = 17.1700; t(i) = 2.7500; i = i + 1 234 y(i) = 10.3900; t(i) = 3.7500; i = i + 1 235 y(i) = 28.4000; t(i) = 1.7500; i = i + 1 236 y(i) = 28.6900; t(i) = 1.7500; i = i + 1 237 y(i) = 81.3000; t(i) = .5000; i = i + 1 238 y(i) = 60.9000; t(i) = .7500; i = i + 1 239 y(i) = 16.6500; t(i) = 2.7500; i = i + 1 240 y(i) = 10.0500; t(i) = 3.7500; i = i + 1 241 y(i) = 28.9000; t(i) = 1.7500; i = i + 1 242 y(i) = 28.9500; t(i) = 1.7500; i = i + 1 243 244 end 245 246 subroutine TaskWorker(ierr) 247 248 PetscErrorCode ierr 249 PetscReal x(n), f(1) 250 PetscMPIInt tag 251 PetscInt index 252 PetscMPIInt status(MPI_STATUS_SIZE) 253 254 tag = IDLE_TAG 255 f = 0.0 256 ! Send check-in message to rank-0 257 PetscCallMPI(MPI_Send(f, one, MPIU_SCALAR, zero, IDLE_TAG, PETSC_COMM_WORLD, ierr)) 258 do while (tag /= DIE_TAG) 259 PetscCallMPI(MPI_Recv(x, nn, MPIU_SCALAR, zero, MPI_ANY_TAG, PETSC_COMM_WORLD, status, ierr)) 260 tag = status(MPI_TAG) 261 if (tag == IDLE_TAG) then 262 PetscCallMPI(MPI_Send(f, one, MPIU_SCALAR, zero, IDLE_TAG, PETSC_COMM_WORLD, ierr)) 263 else if (tag /= DIE_TAG) then 264 index = tag 265 ! Compute local part of residual 266 PetscCall(RunSimulation(x, index, f(1), ierr)) 267 268 ! Return residual to rank-0 269 PetscCallMPI(MPI_Send(f, one, MPIU_SCALAR, zero, tag, PETSC_COMM_WORLD, ierr)) 270 end if 271 end do 272 ierr = 0 273 end 274 275 subroutine RunSimulation(x, i, f, ierr) 276 277 PetscReal x(n), f 278 PetscInt i 279 PetscErrorCode ierr 280 f = y(i) - exp(-x(1)*t(i))/(x(2) + x(3)*t(i)) 281 ierr = 0 282 end 283 284 subroutine StopWorkers(ierr) 285 286 integer checkedin 287 PetscMPIInt status(MPI_STATUS_SIZE) 288 PetscMPIInt source 289 PetscReal f(1), x(n) 290 PetscErrorCode ierr 291 PetscInt i 292 293 checkedin = 0 294 do while (checkedin < size - 1) 295 PetscCallMPI(MPI_Recv(f, one, MPIU_SCALAR, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, status, ierr)) 296 checkedin = checkedin + 1 297 source = status(MPI_SOURCE) 298 do i = 1, n 299 x(i) = 0.0 300 end do 301 PetscCallMPI(MPI_Send(x, nn, MPIU_SCALAR, source, DIE_TAG, PETSC_COMM_WORLD, ierr)) 302 end do 303 ierr = 0 304 end 305 306! -------------------------------------------------------------------- 307! FormFunction - Evaluates the function f(X) and gradient G(X) 308! 309! Input Parameters: 310! tao - the Tao context 311! X - input vector 312! dummy - not used 313! 314! Output Parameters: 315! f - function vector 316 317 subroutine FormFunction(ta, x, f, dummy, ierr) 318 319 Tao ta 320 Vec x, f 321 PetscErrorCode ierr 322 323 PetscInt i, checkedin 324 PetscInt finished_tasks 325 PetscMPIInt next_task 326 PetscMPIInt status(MPI_STATUS_SIZE), tag, source 327 PetscInt dummy 328 329 PetscReal, pointer :: f_v(:), x_v(:) 330 PetscReal fval(1) 331 332 ierr = 0 333 334! Get pointers to vector data 335 PetscCall(VecGetArrayRead(x, x_v, ierr)) 336 PetscCall(VecGetArray(f, f_v, ierr)) 337 338! Compute F(X) 339 if (size == 1) then 340 ! Single processor 341 do i = 1, m 342 PetscCall(RunSimulation(x_v, i, f_v(i), ierr)) 343 end do 344 else 345 ! Multiprocessor main 346 next_task = zero 347 finished_tasks = 0 348 checkedin = 0 349 350 do while (finished_tasks < m .or. checkedin < size - 1) 351 PetscCallMPI(MPI_Recv(fval, one, MPIU_SCALAR, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, status, ierr)) 352 tag = status(MPI_TAG) 353 source = status(MPI_SOURCE) 354 if (tag == IDLE_TAG) then 355 checkedin = checkedin + 1 356 else 357 f_v(tag + 1) = fval(1) 358 finished_tasks = finished_tasks + 1 359 end if 360 if (next_task < m) then 361 ! Send task to worker 362 PetscCallMPI(MPI_Send(x_v, nn, MPIU_SCALAR, source, next_task, PETSC_COMM_WORLD, ierr)) 363 next_task = next_task + one 364 else 365 ! Send idle message to worker 366 PetscCallMPI(MPI_Send(x_v, nn, MPIU_SCALAR, source, IDLE_TAG, PETSC_COMM_WORLD, ierr)) 367 end if 368 end do 369 end if 370 371! Restore vectors 372 PetscCall(VecRestoreArrayRead(x, x_v, ierr)) 373 PetscCall(VecRestoreArray(F, f_v, ierr)) 374 end 375 376 subroutine FormStartingPoint(x) 377 378 Vec x 379 PetscReal, pointer :: x_v(:) 380 PetscErrorCode ierr 381 382 PetscCall(VecGetArray(x, x_v, ierr)) 383 x_v(1) = 0.15 384 x_v(2) = 0.008 385 x_v(3) = 0.01 386 PetscCall(VecRestoreArray(x, x_v, ierr)) 387 end 388end module chwirut2fmodule 389 390program main 391 use chwirut2fmodule 392 implicit none 393! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 394! Variable declarations 395! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 396! 397! See additional variable declarations in the file chwirut2f.h 398 399 PetscErrorCode ierr ! used to check for functions returning nonzeros 400 Vec x ! solution vector 401 Vec f ! vector of functions 402 Tao ta ! Tao context 403 404! Initialize TAO and PETSc 405 PetscCallA(PetscInitialize(ierr)) 406 PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr)) 407 PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) 408 409! Initialize problem parameters 410 call InitializeData() 411 412 if (rank == 0) then 413! Allocate vectors for the solution and gradient 414 PetscCallA(VecCreateSeq(PETSC_COMM_SELF, n, x, ierr)) 415 PetscCallA(VecCreateSeq(PETSC_COMM_SELF, m, f, ierr)) 416 417! The TAO code begins here 418 419! Create TAO solver 420 PetscCallA(TaoCreate(PETSC_COMM_SELF, ta, ierr)) 421 PetscCallA(TaoSetType(ta, TAOPOUNDERS, ierr)) 422 423! Set routines for function, gradient, and hessian evaluation 424 PetscCallA(TaoSetResidualRoutine(ta, f, FormFunction, 0, ierr)) 425 426! Optional: Set initial guess 427 call FormStartingPoint(x) 428 PetscCallA(TaoSetSolution(ta, x, ierr)) 429 430! Check for TAO command line options 431 PetscCallA(TaoSetFromOptions(ta, ierr)) 432! SOLVE THE APPLICATION 433 PetscCallA(TaoSolve(ta, ierr)) 434 435! Free TAO data structures 436 PetscCallA(TaoDestroy(ta, ierr)) 437 438! Free PETSc data structures 439 PetscCallA(VecDestroy(x, ierr)) 440 PetscCallA(VecDestroy(f, ierr)) 441 PetscCallA(StopWorkers(ierr)) 442 443 else 444 PetscCallA(TaskWorker(ierr)) 445 end if 446 447 PetscCallA(PetscFinalize(ierr)) 448end 449!/*TEST 450! 451! build: 452! requires: !complex 453! 454! test: 455! nsize: 3 456! args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5 457! requires: !single 458! 459! 460!TEST*/ 461