[petsc-users] NASM FormFunctionLocal

Rongliang Chen rongliang.chan at gmail.com
Sun May 31 09:05:11 CDT 2015


Hi Barry,

Thanks for your reply.

I have no idea how to pass a vector into the FormFunctionLocal() right 
now. I will try to read the NASM source code and hope can find a way to 
do it.

Best regards,
Rongliang

On 05/31/2015 10:41 AM, Barry Smith wrote:
>    Rongliang,
>
>      It appears the NASM solver creates its own subdomains for which it calls the local function on each one:
>
> PetscErrorCode SNESSetUp_NASM(SNES snes)
> ....
>   if (!nasm->subsnes) {
>      ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
>      if (dm) {
>        nasm->usesdm = PETSC_TRUE;
>        ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
>        if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
>        ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
>
>        ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
>        ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr);
>        for (i=0; i<nasm->n; i++) {
>          ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
>
> these subdomains do not correspond to the "normal" one subdomain per process decomposition of the domain that you expect when you call
>
> ierr = DMDAVecGetArray(info->da,user->temp, &z);CHKERRQ(ierr); // added by rlchen
>
> inside PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
>
> In other words you are expecting that the info->da passed into FormFunctionLocal() is the same as the DA in your original main() program. But with NASM solver it simply is not, it is some other DA that is created by the NASM solver.
>
> If you want to pass some vector like your .temp vector into your local function you will need to understand exactly how the NASM solver works so you can pass the correct vectors in. Simple passing a vector associated with the original DA of the main program won't work.
>
>    Barry
>
>
>
>
>> On May 30, 2015, at 8:02 PM, Rongliang Chen <rongliang.chan at gmail.com> wrote:
>>
>> Hi Barry and Matt,
>>
>> Many thanks for your reply. I only have one DA and I do not think I changed others when I changed the overlap.
>>
>> This error can be reproduced using the attached code (added some lines in the src/snes/examples/tutorials/ex19.c) run with:
>>
>> runex19_1:
>>         -@${MPIEXEC} -n 2 ./ex19 -da_refine 3 -snes_type nasm -da_overlap 2 -snes_monitor
>>
>> The full error messages are followed:
>> ------------------
>> lid velocity = 0.0016, prandtl # = 1, grashof # = 1
>>   0 SNES Function norm 4.066115181565e-02
>> [1]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
>> [1]PETSC ERROR: Arguments are incompatible
>> [1]PETSC ERROR: Vector local size 1300 is not compatible with DMDA local sizes 1400 1728
>>
>> [1]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
>> [1]PETSC ERROR: Petsc Release Version 3.5.2, unknown
>> [1]PETSC ERROR: ./ex19 on a 64bit-debug named rlchen by rlchen Sun May 31 08:55:16 2015
>> [1]PETSC ERROR: Configure options --download-fblaslapack --download-blacs --download-scalapack --download-metis --download-parmetis --download-exodusii --download-netcdf --download-hdf5 --with-64-bit-indices --with-c2html=0 --with-mpi=1 --with-debugging=1 --with-shared-libraries=0
>> [1]PETSC ERROR: #1 DMDAVecGetArray() line 71 in /home/rlchen/soft/petsc-3.5.2/src/dm/impls/da/dagetarray.c
>> [1]PETSC ERROR: #2 FormFunctionLocal() line 275 in /home/rlchen/soft/petsc-3.5.2/src/snes/examples/tutorials/ex19.c
>> [1]PETSC ERROR: #3 SNESComputeFunction_DMDA() line 90 in /home/rlchen/soft/petsc-3.5.2/src/snes/utils/dmdasnes.c
>> [1]PETSC ERROR: #4 SNESComputeFunction() line 2033 in /home/rlchen/soft/petsc-3.5.2/src/snes/interface/snes.c
>> [1]PETSC ERROR: #5 SNESSolve_NEWTONLS() line 174 in /home/rlchen/soft/petsc-3.5.2/src/snes/impls/ls/ls.c
>> [1]PETSC ERROR: #6 SNESSolve() line 3743 in /home/rlchen/soft/petsc-3.5.2/src/snes/interface/snes.c
>> [1]PETSC ERROR: #7 SNESNASMSolveLocal_Private() line 716 in /home/rlchen/soft/petsc-3.5.2/src/snes/impls/nasm/nasm.c
>> [1]PETSC ERROR: #8 SNESSolve_NASM() line 859 in /home/rlchen/soft/petsc-3.5.2/src/snes/impls/nasm/nasm.c
>> [1]PETSC ERROR: #9 SNESSolve() line 3743 in /home/rlchen/soft/petsc-3.5.2/src/snes/interface/snes.c
>> [1]PETSC ERROR: #10 main() line 162 in /home/rlchen/soft/petsc-3.5.2/src/snes/examples/tutorials/ex19.c
>> [1]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
>> application called MPI_Abort(MPI_COMM_WORLD, 75) - process 1
>> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
>> [0]PETSC ERROR: Arguments are incompatible
>> [0]PETSC ERROR: Vector local size 1400 is not compatible with DMDA local sizes 1500 1836
>>
>> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
>> [0]PETSC ERROR: Petsc Release Version 3.5.2, unknown
>> [0]PETSC ERROR: ./ex19 on a 64bit-debug named rlchen by rlchen Sun May 31 08:55:16 2015
>> [0]PETSC ERROR: Configure options --download-fblaslapack --download-blacs --download-scalapack --download-metis --download-parmetis --download-exodusii --download-netcdf --download-hdf5 --with-64-bit-indices --with-c2html=0 --with-mpi=1 --with-debugging=1 --with-shared-libraries=0
>> [0]PETSC ERROR: #1 DMDAVecGetArray() line 71 in /home/rlchen/soft/petsc-3.5.2/src/dm/impls/da/dagetarray.c
>> [0]PETSC ERROR: #2 FormFunctionLocal() line 275 in /home/rlchen/soft/petsc-3.5.2/src/snes/examples/tutorials/ex19.c
>> [0]PETSC ERROR: #3 SNESComputeFunction_DMDA() line 90 in /home/rlchen/soft/petsc-3.5.2/src/snes/utils/dmdasnes.c
>> [0]PETSC ERROR: #4 SNESComputeFunction() line 2033 in /home/rlchen/soft/petsc-3.5.2/src/snes/interface/snes.c
>> [0]PETSC ERROR: #5 SNESSolve_NEWTONLS() line 174 in /home/rlchen/soft/petsc-3.5.2/src/snes/impls/ls/ls.c
>> [0]PETSC ERROR: #6 SNESSolve() line 3743 in /home/rlchen/soft/petsc-3.5.2/src/snes/interface/snes.c
>> [0]PETSC ERROR: #7 SNESNASMSolveLocal_Private() line 716 in /home/rlchen/soft/petsc-3.5.2/src/snes/impls/nasm/nasm.c
>> [0]PETSC ERROR: #8 SNESSolve_NASM() line 859 in /home/rlchen/soft/petsc-3.5.2/src/snes/impls/nasm/nasm.c
>> [0]PETSC ERROR: #9 SNESSolve() line 3743 in /home/rlchen/soft/petsc-3.5.2/src/snes/interface/snes.c
>> [0]PETSC ERROR: #10 main() line 162 in /home/rlchen/soft/petsc-3.5.2/src/snes/examples/tutorials/ex19.c
>> [0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
>> application called MPI_Abort(MPI_COMM_WORLD, 75) - process 0
>>
>> =====================================================================================
>> =   BAD TERMINATION OF ONE OF YOUR APPLICATION PROCESSES
>> =   EXIT CODE: 19200
>> =   CLEANING UP REMAINING PROCESSES
>> =   YOU CAN IGNORE THE BELOW CLEANUP MESSAGES
>> =====================================================================================
>> make: [runex19_1] Error 75 (ignored)
>> -----------------------
>>
>> Best regards,
>> Rongliang
>>
>> On 05/31/2015 01:09 AM, Barry Smith wrote:
>>>> On May 30, 2015, at 5:05 AM, Rongliang Chen <rongliang.chan at gmail.com> wrote:
>>>>
>>>> Hi there,
>>>>
>>>> I am trying to use NASM to solve a nonlinear problem. In my FormFunctionLocal, I need to call DMDAVecGetArray(da, myvec, z). For the "-da_overlap 0" case, the code works fine.
>>>>
>>>> For the "-da_overlap 1" case, when I use the vector "myvec" created by DMCreateLocalVector, it works fine. But if I use the vector "myvec" created by DMCreateGlobalVector, I got the error messages:
>>>> -----------------------
>>>> [6]Arguments are incompatible
>>>> [6]PETSC ERROR: Vector local size 597870 is not compatible with DMDA local sizes 619528 664392
>>> I think you are changing something else at the same time you changed the overlap.  The number ^^^^  619528 is the local size of the global vector and it should remain the same when you change the overlap but below you have TWO other numbers 633144 and 650256  Do you have more than one DMDA?
>>> Run changing absolutely nothing except the overlap.
>>>> ----------------------
>>>>
>>>> And when I switch to "-da_overlap 2", both of the global and local version of "mrvec" got error messages:
>>>> ----------------------------------
>>>> [10]PETSC ERROR: Arguments are incompatible
>>>> [10]PETSC ERROR: Vector local size 589680 is not compatible with DMDA local sizes 633144 678680
>>>>
>>>> [3]PETSC ERROR: Arguments are incompatible
>>>> [3]PETSC ERROR: Vector local size 619528 is not compatible with DMDA local sizes 650256 696540
>>>> -----------------------------------
>>>>
>>>> Can any one tell me how to use DMDAVecGetArray in FormFunctionLocal for NASM? Thanks.
>>>>
>>>> Best regards,
>>>> Rongliang
>>>>
>> <ex19.c>



More information about the petsc-users mailing list